API Reference

Types

SpectroscopyTools.SpectroscopyToolsModule

SpectroscopyTools.jl — General-purpose spectroscopy analysis toolkit.

Provides data types, fitting routines, baseline correction, peak detection, unit conversions, and utility functions for spectroscopic data analysis.

Extracted from QPS.jl to serve as a standalone, reusable foundation.

source
SpectroscopyTools.AbstractSpectroscopyDataType
AbstractSpectroscopyData

Abstract base type for all spectroscopy data types.

Subtypes must implement the following interface:

  • xdata(d) — Primary x-axis data (Vector{Float64})
  • ydata(d) — Signal data (Vector{Float64}) or secondary axis for 2D
  • zdata(d) — Matrix data for 2D types, nothing for 1D
  • xlabel(d) — X-axis label string
  • ylabel(d) — Y-axis or signal label string
  • is_matrix(d) — Whether data is 2D (returns Bool)

Optional (have defaults):

  • source_file(d) — Source filename (default: "")
  • npoints(d) — Number of data points (default: length(xdata(d)), tuple for 2D)
  • title(d) — Display title (default: source_file(d))

This enables uniform handling in data viewers and plotting functions while maintaining semantic field names in each concrete type.

Example

# Works for any spectroscopy data type
function plot_data(data::AbstractSpectroscopyData)
    if is_matrix(data)
        heatmap(xdata(data), ydata(data), zdata(data)')
    else
        lines(xdata(data), ydata(data))
    end
end
source
SpectroscopyTools.TATraceType
TATrace <: AbstractSpectroscopyData

Single-wavelength transient absorption kinetic trace.

Fields

  • time::Vector{Float64}: Time axis (ps)
  • signal::Vector{Float64}: ΔA signal
  • wavelength::Float64: Probe wavelength, NaN if unknown
  • metadata::Dict{Symbol,Any}: Additional info
source
SpectroscopyTools.TASpectrumType
TASpectrum <: AbstractSpectroscopyData

Transient absorption spectrum at a fixed time delay.

Fields

  • wavenumber::Vector{Float64}: Wavenumber axis (cm⁻¹)
  • signal::Vector{Float64}: ΔA signal
  • time_delay::Float64: Time delay (ps), NaN if unknown
  • metadata::Dict{Symbol,Any}: Additional info
source
SpectroscopyTools.TAMatrixType
TAMatrix <: AbstractSpectroscopyData

Two-dimensional transient absorption data (time x wavelength).

Fields

  • time::Vector{Float64}: Time axis (ps)
  • wavelength::Vector{Float64}: Wavelength axis (nm) or wavenumber (cm⁻¹)
  • data::Matrix{Float64}: ΔA signal matrix, size (ntime, nwavelength)
  • metadata::Dict{Symbol,Any}: Additional info

Indexing

matrix[λ=800]     # Extract TATrace at λ ≈ 800 nm
matrix[t=1.0]     # Extract TASpectrum at t ≈ 1.0 ps
source
SpectroscopyTools.TASpectrumFitType
TASpectrumFit

Result of TA spectrum fitting with N peaks of arbitrary lineshape.

Supports any combination of ESA, GSB, and SE peaks, each with its own lineshape model (Gaussian, Lorentzian, pseudo-Voigt).

Access peaks

  • result.peaks — Vector of TAPeak
  • result[i] — i-th peak
  • result[:esa] — first peak with label :esa
  • anharmonicity(result) — GSB center minus ESA center (NaN if not applicable)

Fields

  • peaks::Vector{TAPeak} — Fitted peak parameters
  • offset, rsquared, residuals — Fit metadata
source
SpectroscopyTools.TAPeakType
TAPeak

Information about a single peak in a TA spectrum fit.

Fields

  • label::Symbol — Peak type (:esa, :gsb, :se, :positive, :negative)
  • model::String — Lineshape model name ("gaussian", "lorentzian", "pseudo_voigt")
  • center::Float64 — Peak center position
  • width::Float64 — Width parameter (sigma for gaussian/voigt, fwhm for lorentzian)
  • amplitude::Float64 — Peak amplitude (always positive; sign determined by label)
source

Chirp Correction

SpectroscopyTools.ChirpCalibrationType
ChirpCalibration

Stores the result of chirp detection: detected chirp points, polynomial fit, and detection parameters for reproducibility.

Fields

  • wavelength: Wavelength points where chirp was detected (nm)
  • time_offset: Detected chirp time at each wavelength (ps)
  • poly_coeffs: Polynomial fit coefficients (constant term first, ascending order)
  • poly_order: Polynomial order used
  • reference_λ: Reference wavelength where chirp = 0 (nm)
  • r_squared: Polynomial fit quality
  • metadata: Detection parameters for reproducibility
source
SpectroscopyTools.detect_chirpFunction
detect_chirp(matrix::TAMatrix; kwargs...) -> ChirpCalibration

Detect chirp (GVD) in a broadband TA matrix via cross-correlation (:xcorr) or threshold crossing (:threshold). Returns a ChirpCalibration with polynomial fit.

source
SpectroscopyTools.correct_chirpFunction
correct_chirp(matrix::TAMatrix, cal::ChirpCalibration) -> TAMatrix

Apply chirp correction via cubic spline interpolation, shifting each wavelength column by t_shift(λ) from the calibration polynomial.

source
SpectroscopyTools.subtract_backgroundFunction
subtract_background(matrix::TAMatrix; t_range=nothing) -> TAMatrix

Subtract pre-pump background from a TA matrix by averaging and removing the signal in the baseline region (before pump arrival).

source
subtract_background(m::PLMap; positions=nothing, margin=5) -> PLMap

Subtract a background spectrum from every grid point.

The background is the average CCD spectrum over the reference positions. After subtraction, the intensity map is recomputed from the corrected spectra using the same pixel_range as the original load (if any).

Arguments

  • positions: Vector of (x, y) spatial coordinate tuples (μm) for background reference points. These should be off-flake positions with no PL signal.
  • margin: Number of grid points from each edge used for auto-detection when positions is not given. Auto mode averages the corners of the bottom half of the map (avoids top-row artifacts). Default: 5.

Example

m = load_pl_map("data.lvm"; nx=51, ny=51, step_size=2.16, pixel_range=(950, 1100))

# Explicit background positions
m_bg = subtract_background(m; positions=[(-40, -40), (40, -40), (-40, -20)])

# Auto mode (bottom corners)
m_bg = subtract_background(m)
source
SpectroscopyTools.polynomialFunction
polynomial(cal::ChirpCalibration) -> Function

Return a callable polynomial t_shift = poly(λ) from the calibration. Coefficients are in ascending order: c[1] + c[2]*λ + c[3]*λ² + ...

source

SVD Filtering

SpectroscopyTools.svd_filterFunction
svd_filter(matrix::TAMatrix; n_components::Int=5) -> TAMatrix

Denoise a TA matrix by keeping only the first n_components singular value components. Higher-order components (dominated by noise) are discarded.

This is a standard preprocessing step for broadband TA data. Typical usage: denoise first, then subtract background, detect chirp, and correct chirp.

Arguments

  • matrix::TAMatrix: Input time × wavelength ΔA matrix

Keywords

  • n_components::Int=5: Number of singular value components to retain. Use singular_values to inspect the spectrum and choose.

Returns

A new TAMatrix with filtered data. Metadata includes :svd_filtered => true and :svd_n_components => n_components.

Examples

sv = singular_values(matrix)  # inspect singular value spectrum
filtered = svd_filter(matrix; n_components=3)
source
svd_filter(x::AbstractVector, y::AbstractVector, data::AbstractMatrix;
           n_components::Int=5) -> Matrix{Float64}

Denoise a raw data matrix by keeping only the first n_components singular value components. Returns the filtered matrix.

Arguments

  • x: First axis (e.g., time)
  • y: Second axis (e.g., wavelength)
  • data: Matrix of size (length(x), length(y))

Keywords

  • n_components::Int=5: Number of components to retain
source
SpectroscopyTools.singular_valuesFunction
singular_values(matrix::TAMatrix) -> Vector{Float64}

Return the singular values of the TA data matrix. Inspect these to choose n_components for svd_filter — look for a gap between signal and noise components.

Examples

sv = singular_values(matrix)
# Plot sv to find the elbow, then filter:
filtered = svd_filter(matrix; n_components=3)
source
singular_values(data::AbstractMatrix) -> Vector{Float64}

Return the singular values of a raw data matrix.

source

Cosmic Ray Detection

SpectroscopyTools.CosmicRayResultType
CosmicRayResult

Result of cosmic ray detection on a 1D signal.

Fields

  • indices::Vector{Int} — flagged channel indices
  • count::Int — number of flagged channels
source
SpectroscopyTools.CosmicRayMapResultType
CosmicRayMapResult

Result of cosmic ray detection on a PLMap.

Fields

  • mask::BitArray{3}(nx, ny, n_pixel), true = cosmic ray
  • count::Int — total flagged voxels
  • affected_spectra::Int — number of spectra with at least one cosmic ray
  • channel_counts::Vector{Int} — cosmic ray count per spectral channel
source
SpectroscopyTools.detect_cosmic_raysFunction
detect_cosmic_rays(signal::AbstractVector; threshold=5.0) -> CosmicRayResult

Detect cosmic ray spikes in a 1D spectrum using modified z-scores on first differences (Whitaker-Hayes method).

Computes first differences Δ[k] = signal[k+1] - signal[k], then flags channels where the modified z-score exceeds threshold. The modified z-score uses the Median Absolute Deviation (MAD) for robust scale estimation.

Arguments

  • signal: 1D spectral signal
  • threshold: z-score cutoff (default 5.0). Lower values detect more spikes but may flag real features.

Returns

A CosmicRayResult with the flagged channel indices.

Example

result = detect_cosmic_rays(spectrum; threshold=5.0)
println("Found $(result.count) cosmic rays at indices $(result.indices)")
source
detect_cosmic_rays(m::PLMap; threshold=5.0, pixel_range=nothing, max_spike_width=7) -> CosmicRayMapResult

Detect cosmic ray spikes across all spectra in a PLMap using Most Similar Neighbor (MSN) comparison.

For each pixel, selects the 4-connected neighbor with the highest Pearson correlation as the reference spectrum. The residual (pixel minus MSN) is then tested for positive outliers using the MAD-based robust scale estimate. Channels where the residual exceeds threshold × σ above the median residual are flagged.

Using the most-similar neighbor avoids false positives at spatial boundaries: real spectral differences cancel because the MSN is on the same side of the boundary. Only true cosmic ray spikes — absent from all neighbors — remain.

Additional filters:

  • Width filter: contiguous flagged runs wider than max_spike_width are unflagged — cosmic rays span 1–5 channels; wider runs are spectral variation.
  • Fraction cap: if more than 5% of channels are still flagged after the width filter, all flags for that pixel are cleared.

A global noise floor (median of per-pixel σ estimates) prevents over-sensitivity in spatially homogeneous regions.

Arguments

  • m: PLMap with 3D spectra array (nx, ny, n_pixel)
  • threshold: outlier cutoff in MAD-scaled units (default 5.0)
  • pixel_range: (start, stop) channel indices to analyze. When set, only this subrange of each spectrum is checked for cosmic rays. Channels outside the range are never flagged. Defaults to the pixel_range in m.metadata, or the full spectrum if unset.
  • max_spike_width: maximum width (in channels) of a contiguous flagged run to keep (default 7). Runs wider than this are unflagged as spectral variation.

Returns

A CosmicRayMapResult with a 3D boolean mask and summary statistics.

Example

cr = detect_cosmic_rays(plmap; threshold=5.0)
println("Found $(cr.count) cosmic rays in $(cr.affected_spectra) spectra")
source
SpectroscopyTools.remove_cosmic_raysFunction
remove_cosmic_rays(signal::AbstractVector, result::CosmicRayResult) -> Vector

Replace flagged cosmic ray channels with linearly interpolated values from the nearest non-flagged neighbors.

Returns a new vector; does not mutate the input.

Example

result = detect_cosmic_rays(spectrum)
cleaned = remove_cosmic_rays(spectrum, result)
source
remove_cosmic_rays(m::PLMap, result::CosmicRayMapResult) -> PLMap

Remove cosmic rays from a PLMap using scaled Most Similar Neighbor (MSN) replacement.

For each affected pixel, finds the 4-connected neighbor with the highest Pearson correlation (the MSN). A scale factor is computed from non-flagged channels to match the pixel's intensity level. Flagged channels are then replaced with scale × MSN[k]. This preserves the pixel's overall intensity while removing only the spike shape, avoiding the dark-patch artifacts that raw neighbor median replacement produces.

Falls back to spectral interpolation (1D linear from nearest non-flagged channels) for edge pixels with no spatial neighbors.

Returns a new PLMap with recomputed intensity.

Example

cr = detect_cosmic_rays(plmap)
cleaned = remove_cosmic_rays(plmap, cr)
source

PL / Raman Mapping

SpectroscopyTools.PLMapType
PLMap <: AbstractSpectroscopyData

Photoluminescence intensity map from a CCD raster scan.

A 2D spatial grid where each point has a full CCD spectrum. The intensity field holds the integrated PL signal at each position; the full spectra are stored in spectra for extraction at individual positions.

Fields

  • intensity::Matrix{Float64} — Integrated PL intensity (nx, ny)
  • spectra::Array{Float64,3} — Raw CCD counts (nx, ny, n_pixel)
  • x::Vector{Float64} — Spatial x positions (μm)
  • y::Vector{Float64} — Spatial y positions (μm)
  • pixel::Vector{Float64} — Pixel indices (or wavelength if calibrated)
  • metadata::Dict{String,Any} — Source file, grid dims, step size, etc.
source
SpectroscopyTools.extract_spectrumFunction
extract_spectrum(m::PLMap, ix::Int, iy::Int) -> NamedTuple

Extract the CCD spectrum at grid index (ix, iy).

Returns (pixel=..., signal=..., x=..., y=...).

source
extract_spectrum(m::PLMap; x::Real, y::Real) -> NamedTuple

Extract the CCD spectrum at the spatial position nearest to (x, y).

Returns (pixel=..., signal=..., x=..., y=..., ix=..., iy=...).

source
SpectroscopyTools.integrated_intensityFunction
integrated_intensity(m::PLMap; pixel_range=nothing) -> Matrix{Float64}

Compute the integrated intensity at each grid point.

Without pixel_range, returns m.intensity (the precomputed full-range sum). With pixel_range, sums m.spectra[:, :, p1:p2] over the given pixel window.

Arguments

  • pixel_range: (start, stop) pixel indices to integrate over. Falls back to the pixel_range stored in metadata, or uses m.intensity if unset.
source
SpectroscopyTools.intensity_maskFunction
intensity_mask(m::PLMap; pixel_range=nothing, threshold=0.05, exclude=nothing) -> NamedTuple

Compute a boolean mask over the PLMap grid based on integrated intensity.

Grid points where the integrated intensity is below threshold fraction of the intensity range are excluded (false). Spatial regions listed in exclude are also set to false regardless of intensity. This is useful for filtering out off-sample regions and known artifacts (e.g., substrate bands) before batch operations like fitting.

The threshold is computed as: cutoff = min + threshold × (max - min).

Arguments

  • pixel_range: (start, stop) pixel indices for computing integrated intensity. If nothing, uses the full precomputed m.intensity.
  • threshold: Fraction of the intensity range (0.0 to 1.0). Default 0.05 (5%).
  • exclude: Spatial regions to exclude. Each region is a tuple of x and y ranges in spatial coordinates (μm): ((x_min, x_max), (y_min, y_max)). Accepts a single region or a vector of regions. Grid points whose spatial coordinates fall within any excluded region are masked out.

Returns

A named tuple with fields:

  • mask::BitMatrixtrue for included grid points
  • n_included::Int — number of included points
  • n_total::Int — total number of grid points
  • intensity_min::Float64 — minimum integrated intensity
  • intensity_max::Float64 — maximum integrated intensity
  • cutoff::Float64 — absolute intensity cutoff value

Examples

m = load_pl_map("scan.lvm"; nx=51, ny=51)
m = subtract_background(m)

# Threshold only
result = intensity_mask(m; pixel_range=(950, 1100), threshold=0.1)

# Exclude a substrate band at the top of the map
result = intensity_mask(m; threshold=0.05,
    exclude=((-Inf, Inf), (40.0, Inf)))

# Multiple exclusion regions
result = intensity_mask(m; threshold=0.05,
    exclude=[
        ((-Inf, Inf), (40.0, Inf)),     # top substrate band
        ((-50.0, -40.0), (-50.0, -40.0))  # noisy corner
    ])
source
SpectroscopyTools.peak_centersFunction
peak_centers(m::PLMap; pixel_range=nothing, threshold=0.05) -> Matrix{Float64}

Compute the centroid (intensity-weighted average pixel) at each grid point.

Returns a (nx, ny) matrix of peak center positions in pixel units. Grid points where the PL intensity (m.intensity) is below threshold × the map maximum are set to NaN (renders as transparent in heatmaps with nan_color=:transparent).

Masking uses the PLMap's intensity field — the integrated PL signal already stored in the map. This produces clean masks that match the intensity heatmap: off-flake regions (low PL signal) are transparent, on-flake regions show centroid positions.

Arguments

  • pixel_range: (start, stop) pixel range to compute the centroid over. Falls back to the pixel_range stored in metadata, or all pixels if unset.
  • threshold: Fraction of the maximum PL intensity below which a grid point is masked as NaN. Default 0.05 (5%). Set to 0 to disable masking.

Example

m = load_pl_map("scan.lvm"; nx=51, ny=51, pixel_range=(950, 1100))
m = subtract_background(m)
centers = peak_centers(m)
heatmap(m.x, m.y, centers; colormap=:viridis, nan_color=:transparent)
source

Peak Fitting

SpectroscopyTools.fit_peaksFunction
fit_peaks(x, y; kwargs...) -> MultiPeakFitResult
fit_peaks(spec::AbstractSpectroscopyData, region; kwargs...) -> MultiPeakFitResult
fit_peaks(spec::AbstractSpectroscopyData; kwargs...) -> MultiPeakFitResult

Fit one or more peaks in spectroscopy data.

source
SpectroscopyTools.fit_ta_spectrumFunction
fit_ta_spectrum(spec::TASpectrum; kwargs...) -> TASpectrumFit

Fit a transient absorption spectrum with N peaks of arbitrary lineshape.

Uses find_peaks to automatically detect initial peak positions from the data, so multiple well-separated peaks of the same type (e.g., three GSB peaks for W(CO)₆) are initialized correctly.

Keywords

  • peaks=[:esa, :gsb] — Peak types. Each element is either a Symbol (:esa, :gsb, :se, :positive, :negative) or a (Symbol, Function) tuple specifying label and lineshape model.
  • model=gaussian — Default lineshape for peaks specified as symbols only.
  • region=nothing — Optional (x_min, x_max) fitting region.
  • fit_offset=false — Whether to fit a constant offset.
  • p0=nothing — Manual initial parameter vector. Overrides automatic detection.

Peak signs

  • :esa, :positive → +1 (positive ΔA)
  • :gsb, :se, :negative → -1 (negative ΔA)

Examples

# Default: 1 Gaussian ESA + 1 Gaussian GSB
result = fit_ta_spectrum(spec)

# Three GSB peaks (e.g., W(CO)₆ carbonyl stretches)
result = fit_ta_spectrum(spec; peaks=[:esa, :esa, :esa, :gsb, :gsb, :gsb])

# Per-peak lineshapes
result = fit_ta_spectrum(spec; peaks=[(:esa, lorentzian), (:gsb, gaussian)])

# Access results
result[:esa].center      # first ESA peak
result[2].center         # second peak by index
anharmonicity(result)    # GSB - ESA center (only if exactly 1 of each)
predict(result, ν)       # full fitted curve
predict_peak(result, 1)  # single peak contribution
source
SpectroscopyTools.MultiPeakFitResultType
MultiPeakFitResult

Result of multi-peak fitting (1 to N peaks with polynomial baseline).

Supports indexing by peak number: result[1] returns first peak's PeakFitResult.

source
SpectroscopyTools.anharmonicityFunction
anharmonicity(fit::TASpectrumFit) -> Float64

Compute the anharmonicity (GSB center - ESA center) from a TA spectrum fit. Returns NaN if there is not exactly one ESA and one GSB peak.

source

Peak Detection

SpectroscopyTools.PeakInfoType
PeakInfo

Information about a detected peak.

Fields

  • position::Float64 — Peak center in x-units
  • intensity::Float64 — Peak height
  • prominence::Float64 — How much the peak stands out from surrounding baseline
  • width::Float64 — Full width at half prominence (FWHP) in x-units
  • bounds::Tuple{Float64,Float64} — Left and right edges at half prominence
  • index::Int — Index in the original data array
source
SpectroscopyTools.find_peaksFunction
find_peaks(y; kwargs...) -> Vector{PeakInfo}
find_peaks(x, y; kwargs...) -> Vector{PeakInfo}

Find peaks in spectroscopic data.

Keyword Arguments

  • min_prominence::Real=0.05 — Minimum prominence as fraction of data range
  • min_width::Real=0 — Minimum peak width in x-units
  • max_width::Real=Inf — Maximum peak width in x-units
  • min_height::Real=-Inf — Minimum peak height (absolute)
  • window::Int=1 — Comparison window for local maxima detection
  • mode::Symbol=:maxima:maxima for peaks, :minima for valleys (e.g. 2nd derivative)
  • baseline::Union{Symbol,Nothing}=nothing — Apply baseline correction before peak detection
  • baseline_kw::NamedTuple=NamedTuple() — Keyword arguments for baseline correction
source

Exponential Decay Fitting

SpectroscopyTools.fit_exp_decayFunction
fit_exp_decay(trace::TATrace; n_exp=1, irf=false, irf_width=0.15, t_start=0.0, t_range=nothing)

Fit exponential decay to a transient absorption trace.

Arguments

  • trace: TATrace
  • n_exp: Number of exponential components (default 1)
  • irf: Include IRF convolution (default false)
  • irf_width: Initial guess for IRF σ in ps (default 0.15)
  • t_start: Start time for fitting when irf=false (default 0.0)
  • t_range: Optional (tmin, tmax) to restrict fit region

Returns

  • n_exp=1: ExpDecayFit
  • n_exp>1: MultiexpDecayFit
source
SpectroscopyTools.fit_globalFunction
fit_global(traces::Vector{TATrace}; n_exp=1, irf_width=0.15, labels=nothing) -> GlobalFitResult

Fit multiple traces simultaneously with shared time constant(s) τ.

Supports single-exponential (n_exp=1) and multi-exponential (n_exp>1) global analysis. All traces share the same time constants, IRF width, and time zero, while amplitudes and offsets are fitted per-trace.

Keywords

  • n_exp::Int=1: Number of exponential components
  • irf_width::Float64=0.15: Initial guess for IRF σ in ps
  • labels=nothing: Optional trace labels (defaults to "Trace 1", "Trace 2", ...)

Returns

GlobalFitResult with shared taus and per-trace amplitudes matrix.

Examples

# Single-exponential global fit
result = fit_global([trace_esa, trace_gsb])

# Multi-exponential with 2 shared time constants
result = fit_global([trace1, trace2, trace3]; n_exp=2)
report(result)
source
fit_global(matrix::TAMatrix; n_exp=1, irf_width=0.15, λ=nothing) -> GlobalFitResult

Global analysis of a TAMatrix, extracting traces at each wavelength.

Returns a GlobalFitResult with the wavelengths field populated, enabling decay-associated spectra (DAS) via das(result).

Keywords

  • n_exp::Int=1: Number of exponential components
  • irf_width::Float64=0.15: Initial guess for IRF σ in ps
  • λ=nothing: Specific wavelengths to fit. If nothing, fits all wavelengths.

Examples

result = fit_global(matrix; n_exp=2)
report(result)

# Get decay-associated spectra
spectra = das(result)  # n_exp × n_wavelengths matrix
source
SpectroscopyTools.ExpDecayFitType
ExpDecayFit

Result of exponential decay fitting with instrument response function convolution.

Fields

  • amplitude, tau, t0, sigma, offset, signal_type, residuals, rsquared
source
SpectroscopyTools.MultiexpDecayFitType
MultiexpDecayFit

Result of multi-exponential decay fitting (n >= 1 components).

Fields

  • taus::Vector{Float64}: Time constants (sorted fast->slow)
  • amplitudes::Vector{Float64}: Corresponding amplitudes
  • t0, sigma, offset, signal_type, residuals, rsquared

Derived properties

  • n_exp(fit): Number of exponential components
  • weights(fit): Relative amplitude weights (normalized to 100%)
source
SpectroscopyTools.GlobalFitResultType
GlobalFitResult

Result of global fitting multiple traces with shared time constants.

Supports single-exponential (nexp=1) and multi-exponential (nexp>1) global analysis with shared τ values and per-trace amplitudes.

Fields

  • taus::Vector{Float64}: Shared time constants (sorted fast→slow)
  • sigma::Float64: Shared IRF width
  • t0::Float64: Shared time zero
  • amplitudes::Matrix{Float64}: Per-trace amplitudes (ntraces × nexp)
  • offsets::Vector{Float64}: Per-trace offsets
  • labels::Vector{String}: Trace labels
  • wavelengths::Union{Nothing, Vector{Float64}}: Wavelength axis (from TAMatrix input)
  • rsquared::Float64: Global R²
  • rsquared_individual::Vector{Float64}: Per-trace R²
  • residuals::Vector{Vector{Float64}}: Per-trace residuals

Derived properties

  • n_exp(fit): Number of exponential components
  • das(fit): Decay-associated spectra (requires TAMatrix input)
source
SpectroscopyTools.dasFunction
das(fit::GlobalFitResult) -> Matrix{Float64}

Return the decay-associated spectra (DAS) as an n_exp × n_wavelengths matrix. Each row is the amplitude spectrum for one time constant.

Requires that the fit was performed on a TAMatrix (wavelengths must be available).

source

Baseline Correction

SpectroscopyTools.rubberband_baselineFunction
rubberband_baseline(x, y)

Rubber band baseline correction using the lower convex hull.

Equivalent to stretching a rubber band under the spectrum – the baseline follows the lower envelope. Used in OPUS for FTIR baseline correction.

Returns the baseline y-values (same length as input).

source
SpectroscopyTools.imodpoly_baselineFunction
imodpoly_baseline(x, y; poly_order=4, maxiter=100, tol=1e-3) -> Vector

Improved Modified Polynomial baseline (Lieber & Mahadevan-Jansen, 2003).

Iteratively fits a polynomial to the spectrum, removing points that are above the fit by more than one standard deviation of the residuals. Converges when the fit changes less than tol between iterations.

source
SpectroscopyTools.rolling_ball_baselineFunction
rolling_ball_baseline(y; half_window=50, smooth_half_window=nothing) -> Vector

Rolling ball baseline estimation using morphological operations.

Applies erosion (rolling minimum) then dilation (rolling maximum) to estimate the baseline envelope, followed by moving-average smoothing.

Arguments

  • half_window::Int=50: Half-window size for morphological operations.
  • smooth_half_window::Int: Half-window for final smoothing (default: half_window).
source
SpectroscopyTools.correct_baselineFunction
correct_baseline(y; method=:arpls, kwargs...) -> NamedTuple

Correct baseline and return both corrected spectrum and baseline.

Returns (y=corrected, baseline=baseline).

source
correct_baseline(x, y; kwargs...) -> NamedTuple

Version that also returns x values for convenience. Methods that use x-values (rubberband) receive real x, not dummy indices.

source

Spectroscopy Utilities

SpectroscopyTools.subtract_spectrumFunction
subtract_spectrum(sample, reference; scale=1.0, interpolate=false)

Subtract a reference spectrum from a sample spectrum.

Accepts AbstractSpectroscopyData types (uses xdata/ydata interface) or any objects with .x and .y fields.

Returns (x=..., y=...) NamedTuple.

source
SpectroscopyTools.npointsFunction
npoints(d::AbstractSpectroscopyData) -> Int
npoints(d::TAMatrix) -> Tuple{Int,Int}

Return the number of data points. For 1D data, returns an Int. For 2D data (TAMatrix), returns (n_time, n_wavelength).

source
SpectroscopyTools.titleFunction
title(d::AbstractSpectroscopyData) -> String

Return a display title for the data. Defaults to source_file(d).

source

Unit Conversions