LegendSpecFits

Also see LegendSpecFits on GitHub and the full package documentation.

Modules

Types and constants

    Functions and macros

    Documentation

    LegendSpecFits._get_model_countsMethod
    _get_model_counts(f_fit::Base.Callable,v_ml::NamedTuple,bin_centers::StepRangeLen,bin_widths::StepRangeLen)

    aux. function to get modelled peakshape based on histogram binning and best-fit parameter

    source
    LegendSpecFits.aoe_compton_background_peakshapeMethod
    aoe_compton_background_peakshape(
        x::Real, μ::Real, σ::Real,
        background::Real, δ::Real
    )

    Describes the background shape of a typical A/E Compton peak in a detector as a step like background for MSE events.

    source
    LegendSpecFits.aoe_compton_peakshapeMethod
    aoe_compton_peakshape(
        x::Real, μ::Real, σ::Real, n::Real,
        background::Real, δ::Real
    )

    Describes the shape of a typical A/E Compton peak in a detector as a gaussian SSE peak and a step like background for MSE events.

    source
    LegendSpecFits.array_to_tupleMethod
    array_to_tuple(a::AbstractArray, as_nt::NamedTuple)

    Return a NamedTuple with the values of a in the order given by fieldnames(as_nt).

    source
    LegendSpecFits.background_peakshapeMethod
    background_peakshape(
        x::Real, μ::Real, σ::Real, n::Real,
        skew_fraction::Real, skew_width::Real,
    )

    Describes the background part of the shape of a typical gamma peak in a detector.

    source
    LegendSpecFits.baseline_qcMethod
    baseline_qc(data, qc_config)

    Perform simple Gaussian fits on the baseline disitrbutions for a given data set.

    Returns

    • result: Namedtuple containing the cut and fit results for each baseline variable
    • report: Namedtuple containing plotable objects.
    source
    LegendSpecFits.chi2fitMethod
    fit_chisq(x::AbstractVector{<:Real},y::AbstractVector{<:Real},yerr::AbstractVector{<:Real}, f_fit::Function;pull_t::Vector{<:NamedTuple} = fill(NamedTuple(), first(methods(f_fit)).nargs - 2), v_init::Vector = [])

    Least square fit with chi2 minimization

    Input:

    • x : x-values
    • y : y-values
    • yerr : 1 sigma uncertainty on y
    • ffit : fit/model function. e.g. for a linear function: flin(x,p1,p2) = p1 .* x .+ p2

    The numer of fit parameter is determined with first(methods(f_fit)).nargs - 2. That's why it's important that ffit has the synthax f(x,arg1,arg2,arg3,...) pullt : pull term, a vector of NamedTuple with fields mean and std. A Gaussian pull term is added to the chi2 function to account for systematic uncertainties. If left blank, no pull term is used. v_init : initial value for fit parameter optimization. If left blank, the initial value is set to 1 or guessed roughly for all fit parameters

    Return:

    • result : NamedTuple with the optimized fit parameter and the fit function
    • report:
    source
    LegendSpecFits.ctc_energyFunction
    ctc_energy(e::Array{T}, qdrift::Array{T}, peak::T, window::T) where T<:Real

    Correct for the drift time dependence of the energy by minimizing the ratio of the FWHM and the peak height of the peak around peak in e with a cut window of window. The drift time dependence is given by qdrift.

    Returns

    * `peak`: peak position
    * `window`: window size
    * `fct`: correction factor
    * `fwhm_before`: FWHM before correction
    * `fwhm_after`: FWHM after correction
    * `func`: function to correct energy
    * `func_generic`: generic function to correct energy
    source
    LegendSpecFits.cut_single_peakMethod
    cut_single_peak(x::Array, min_x::Float64, max_x::Float64,; n_bins::Int=15000, relative_cut::Float64=0.5)

    Cut out a single peak from the array x between min_x and max_x. The number of bins is the number of bins to use for the histogram. The relative cut is the fraction of the maximum counts to use for the cut.

    Returns

    * `max`: maximum position of the peak
    * `low`: lower edge of the cut peak
    * `high`: upper edge of the cut peak
    source
    LegendSpecFits.estimate_fwhmMethod
    estimate_fwhm(v::NamedTuple, v_err::NamedTuple)

    Get the FWHM of a peak from the fit parameters.

    Returns

    * `fwhm`: the FWHM of the peak
    source
    LegendSpecFits.estimate_single_peak_statsMethod
    estimate_single_peak_stats(h::Histogram, calib_type::Symbol=:th228)

    Estimate statistics/parameters for a single peak in the given histogram h.

    h must only contain a single peak. The peak should have a Gaussian-like shape. calib_type specifies the calibration type. Currently only :th228 is implemented. If you want get the peak statistics for a PSD calibration, use :psd.

    Returns

    NamedTuple with the fields * peak_pos: estimated position of the peak (in the middle of the peak) * peak_fwhm: full width at half maximum (FWHM) of the peak * peak_sigma: estimated standard deviation of the peak * peak_counts: estimated number of counts in the peak * mean_background: estimated mean background value

    source
    LegendSpecFits.estimate_single_peak_stats_psdMethod
    estimate_single_peak_stats_psd(h::Histogram{T}) where T<:Real

    Estimate peak parameters for a single peak in a histogram using the maximum, the FWHM and the area of the peak.

    Returns

    * `peak_pos`: Position of the peak
    * `peak_fwhm`: Full width at half maximum of the peak
    * `peak_sigma`: Standard deviation of the peak
    * `peak_counts`: Counts of the peak
    * `mean_background`: Mean background of the peak
    source
    LegendSpecFits.ex_step_gaussMethod
    ex_step_gauss(x::Real, l::Real, k::Real, t::Real, d::Real)

    Evaluates an extended step gauss model at x with parameters l, k, t and d.

    source
    LegendSpecFits.expand_varsMethod
    expand_vars(v::NamedTuple)::StructArray

    Expand all fields in v (scalars or arrays) to same array size and return a StructArray.

    source
    LegendSpecFits.exponential_decayMethod
    exponential_decay(x::Real, amplitude::Real, decay::Real, offset::Real)

    Evaluates an exponential decay function at x with parameters amplitude, decay and offset.

    source
    LegendSpecFits.f_optimize_ctcMethod
    f_optimize_ctc(fct, e, qdrift, bin_width)

    Calculate the ratio of the FWHM and the peak height of the peak around peak in e with a cut window of window. The drift time dependence is given by e_ctc = e + fct * qdrift.

    Returns

    * `fwhm / p_height`: FWHM of the peak divided by peak height
    source
    LegendSpecFits.fit_aoe_comptonMethod
    fit_aoe_compton(peakhists::Array, peakstats::StructArray, compton_bands::Array{T}) where T<:Real

    Fit the A/E Compton bands using the f_aoe_compton function consisting of a gaussian SSE peak and a step like background for MSE events.

    Returns

    * `result`: Dict of NamedTuples of the fit results containing values and errors for each compton band
    * `report`: Dict of NamedTuples of the fit report which can be plotted for each compton band
    source
    LegendSpecFits.fit_aoe_correctionsMethod

    fitaoecorrections(e::Array{<:Unitful.Energy{<:Real}}, μ::Array{<:Real}, σ::Array{<:Real})

    Fit the corrections for the AoE value of the detector.

    Returns

    • e: Energy values
    • μ: Mean values
    • σ: Sigma values
    • μ_scs: Fit result for the mean values
    • f_μ_scs: Fit function for the mean values
    • σ_scs: Fit result for the sigma values
    • f_σ_scs: Fit function for the sigma values
    source
    LegendSpecFits.fit_calibrationMethod
    fit_calibration(n_poly::Int, µ::AbstractVector{<:Union{Real,Measurement{<:Real}}}, peaks::AbstractVector{<:Union{Real,Measurement{<:Real}}}; pull_t::Vector{<:NamedTuple}=fill(NamedTuple(), n_poly+1), v_init::Vector = [], uncertainty::Bool=true )

    Fit the calibration lines with polynomial function of npoly order npoly == 1 -> linear function n_poly == 2 -> quadratic function

    Returns

    * `result`: NamedTuple with the following fields
        * `par`: best-fit parameters
        * `gof`: godness of fit
    * `report`:
    source
    LegendSpecFits.fit_enc_sigmasMethod
    fit_enc_sigmas(enc_grid::Matrix{T}, enc_grid_rt::StepRangeLen{Quantity{<:T}, Base.TwicePrecision{Quantity{<:T}}, Base.TwicePrecision{Quantity{<:T}}, Int64}, min_enc::T, max_enc::T, nbins::Int64, rel_cut_fit::T) where T<:Real

    Fit the ENC values in enc_grid for each RT in enc_grid_rt with a Gaussian and return the optimal RT and the corresponding ENC value.

    Arguments

    • enc_grid: 2D array of ENC values for each RT in enc_grid_rt
    • enc_grid_rt: 1D array of RT values for which the ENC values in enc_grid are calculated
    • min_enc: minimum ENC value to consider for the fit
    • max_enc: maximum ENC value to consider for the fit
    • nbins: number of bins to use for the histogram of ENC values
    • rel_cut_fit: relative cut value to use for the fit

    Returns

    • rt: optimal RT value
    • min_enc: corresponding ENC value
    source
    LegendSpecFits.fit_fwhmMethod
    fitFWHM(fit_fwhm(peaks::Vector{T}, fwhm::Vector{T}) where T<:Real

    Fit the FWHM of the peaks to a quadratic function.

    Returns

    * `qbb`: the FWHM at 2039 keV
    * `err`: the uncertainties of the fit parameters
    * `v`: the fit result parameters
    * `f_fit`: the fitted function
    source
    LegendSpecFits.fit_fwhm_ft_fepMethod

    fitfwhmftfep(egrid::Matrix, egridft::StepRangeLen{Quantity{<:T}, Base.TwicePrecision{Quantity{<:T}}, Base.TwicePrecision{Quantity{<:T}}, Int64}, rt::Unitful.RealOrRealQuantity, mine::T, maxe::T, nbins::Int64, relcutfit::T; default_ft::Quantity{T}=3.0u"µs") where {T <:Real}

    Fit the FWHM values in e_grid for each FT in e_grid_ft with a Gamma Peakshape and return the optimal FT and the corresponding FWHM value. The cut values cut for each flat-top time a window for better histogramming.

    Arguments

    • e_grid: 2D array of energy values for each FT in e_grid_ft
    • e_grid_ft: 1D array of FT values for which the FWHM values in e_grid are calculated
    • rt: RT value for which the FWHM values in e_grid are calculated
    • min_e: minimum energy value to consider for the fit
    • max_e: maximum energy value to consider for the fit
    • nbins: number of bins to use for the histogram of energy values
    • rel_cut_fit: relative cut value to use for the fit

    Returns

    • ft: optimal FT value
    • min_fwhm: corresponding FWHM value
    source
    LegendSpecFits.fit_half_centered_trunc_gaussMethod
    fit_half_centered_trunc_gauss(x::Array, cuts::NamedTuple{(:low, :high, :max), Tuple{Float64, Float64, Float64}})

    Fit a single truncated Gaussian to the data x between cut.low and cut.high. The peak center is fixed at μ and the peak is cut in half either in the left or right half.

    Returns report and result` with:

    * `f_fit`: fitted function
    * `μ`: mean of the Gaussian
    * `μ_err`: error of the mean
    * `σ`: standard deviation of the Gaussian
    * `σ_err`: error of the standard deviation
    * `n`: number of counts in the peak
    source
    LegendSpecFits.fit_half_trunc_gaussMethod
    fit_half_centered_trunc_gauss(x::Array, cuts::NamedTuple{(:low, :high, :max), Tuple{Float64, Float64, Float64}})

    Fit a single truncated Gaussian to the data x between cut.low and cut.high. The peak center is fixed at μ and the peak is cut in half either in the left or right half.

    Returns report and result with:

    * `f_fit`: fitted function
    * `μ`: mean of the Gaussian
    * `μ_err`: error of the mean
    * `σ`: standard deviation of the Gaussian
    * `σ_err`: error of the standard deviation
    * `n`: number of counts in the peak
    source
    LegendSpecFits.fit_peaksMethod
    fit_peaks(peakhists::Array, peakstats::StructArray, th228_lines::Array,; calib_type::Symbol=:th228, uncertainty::Bool=true, low_e_tail::Bool=true)

    Perform a fit of the peakshape to the data in peakhists using the initial values in peakstats to the calibration lines in th228_lines.

    Returns

    * `peak_fit_plots`: array of plots of the peak fits
    * `return_vals`: dictionary of the fit results
    source
    LegendSpecFits.fit_peaks_combinedMethod
    fit_peaks_combined(peakhists::Array, peakstats::StructArray, th228_lines::Array{T},; calib_type::Symbol=:th228, uncertainty::Bool=true, fixed_position::Bool=false) where T<:Real

    Fit the peaks of a histogram to a combined peakshape function while sharring parameters between peaks.

    source
    LegendSpecFits.fit_sg_wlMethod
    fit_sg_wl(dep_sep_data, a_grid_wl_sg, optimization_config)

    Fit the SG window length for the SEP data and return the optimal window length and the corresponding survival fraction.

    Arguments

    • dep_sep_data: NamedTuple with the DEP and SEP data
    • a_grid_wl_sg: range of window lengths to sweep through
    • optimization_config: configuration dictionary

    Returns

    • result: optimal window length and corresponding survival fraction
    • report: report with all window lengths and survival fractions
    source
    LegendSpecFits.fit_single_aoe_comptonMethod
    fit_single_aoe_compton(h::Histogram, ps::NamedTuple{(:peak_pos, :peak_fwhm, :peak_sigma, :peak_counts, :mean_background, :μ, :σ), NTuple{7, T}}; uncertainty::Bool=true) where T<:Real

    Perform a fit of the peakshape to the data in h using the initial values in ps while using the f_aoe_compton function consisting of a gaussian SSE peak and a step like background for MSE events.

    Returns

    * `result`: NamedTuple of the fit results containing values and errors
    * `report`: NamedTuple of the fit report which can be plotted
    source
    LegendSpecFits.fit_single_peak_th228Method
    fit_single_peak_th228(h::Histogram, ps::NamedTuple{(:peak_pos, :peak_fwhm, :peak_sigma, :peak_counts, :mean_background), NTuple{5, T}};, uncertainty::Bool=true, fixed_position::Bool=false, low_e_tail::Bool=true) where T<:Real

    Perform a fit of the peakshape to the data in h using the initial values in ps while using the gamma_peakshape with low-E tail. Also, FWHM is calculated from the fitted peakshape with MC error propagation. The peak position can be fixed to the value in ps by setting fixed_position=true. If the low-E tail should not be fitted, it can be disabled by setting low_e_tail=false.

    Returns

    * `result`: NamedTuple of the fit results containing values and errors
    * `report`: NamedTuple of the fit report which can be plotted
    source
    LegendSpecFits.fit_single_trunc_gaussMethod
    fit_single_trunc_gauss(x::Array, cuts::NamedTuple{(:low, :high, :max), Tuple{Float64, Float64, Float64}})

    Fit a single truncated Gaussian to the data x between min_x and max_x. Returns report and resultwith: *ffit: fitted function *μ: mean of the Gaussian *μerr: error of the mean *σ: standard deviation of the Gaussian *σ_err: error of the standard deviation *n`: number of counts in the peak

    source
    LegendSpecFits.fit_subpeaks_th228Method
    fit_single_peak_th228(h::Histogram, ps::NamedTuple{(:peak_pos, :peak_fwhm, :peak_sigma, :peak_counts, :mean_background), NTuple{5, T}};, uncertainty::Bool=true, fixed_position::Bool=false, low_e_tail::Bool=true) where T<:Real

    Perform a simultaneous fit of two peaks (h_survived and h_cut) that together would form a histogram h, from which the result h_result was already determined using fit_single_peak_th228. Also, FWHM is calculated from the fitted peakshape with MC error propagation. The peak position can be fixed to the value in ps by setting fixed_position=true. If the low-E tail should not be fitted, it can be disabled by setting low_e_tail=false.

    Returns

    * `result`: NamedTuple of the fit results containing values and errors, in particular the signal survival fraction `sf` and the background survival frachtion `bsf`.
    * `report`: NamedTuple of the fit report which can be plotted
    source
    LegendSpecFits.gamma_peakshapeMethod
    gamma_peakshape(
        x::Real, μ::Real, σ::Real, n::Real,
        step_amplitude::Real, skew_fraction::Real, skew_width::Real,
        background::Real
    )

    Describes the shape of a typical gamma peak in a detector with a flat background.

    source
    LegendSpecFits.generate_aoe_compton_bandsMethod
    generate_aoe_compton_bands(aoe::Vector{<:Real}, e::Vector{<:T}, compton_bands::Vector{<:T}, compton_window::T) where T<:Unitful.Energy{<:Real}

    Generate histograms for the A/E Compton bands and estimate peak parameters. The compton bands are cutted out of the A/E spectrum and then binned using the Freedman-Diaconis Rule. For better performance the binning is only done in the area around the peak. The peak parameters are estimated using the estimate_single_peak_stats_psd function.

    Returns

    * `peakhists`: Array of histograms for each compton band
    * `peakstats`: StructArray of peak parameters for each compton band
    * `min_aoe`: Array of minimum A/E values for each compton band
    * `max_aoe`: Array of maximum A/E values for each compton band
    * `mean_peak_pos`: Mean peak position of all compton bands
    * `std_peak_pos`: Standard deviation of the peak position of all compton bands
    * `simple_pars_aoe_μ`: Simple curve fit parameters for the peak position energy depencence
    * `simple_pars_error_aoe_μ`: Simple curve fit parameter errors for the peak position energy depencence
    * `simple_pars_aoe_σ`: Simple curve fit parameters for the peak sigma energy depencence
    * `simple_pars_error_aoe_σ`: Simple curve fit parameter errors for the peak sigma energy depencence
    source
    LegendSpecFits.get_aoe_cutMethod
    get_aoe_cut(aoe::Vector{<:Unitful.RealOrRealQuantity}, e::Vector{<:T},; dep::T=1592.53u"keV", window::Vector{T}=[12.0, 10.0]u"keV", dep_sf::Float64=0.9, cut_search_interval::Tuple{Quantity{<:Real}, Quantity{<:Real}}=(-25.0u"keV^-1", 0.0u"keV^-1"), rtol::Float64=0.001, bin_width_window::T=3.0u"keV", fixed_position::Bool=true, sigma_high_sided::Float64=NaN, uncertainty::Bool=true, maxiters::Int=200) where T<:Unitful.Energy{<:Real}

    Get the AoE cut value for a given dep and window size while performing a peak fit with fixed position. The AoE cut value is determined by finding the cut value for which the number of counts after the cut is equal to dep_sf times the number of counts before the cut. The algorhithm utilizes a root search algorithm to find the cut value with a relative tolerance of rtol.

    Returns

    • cut: AoE cut value
    • n0: Number of counts before the cut
    • nsf: Number of counts after the cut
    source
    LegendSpecFits.get_centered_gaussian_window_cutMethod
    get_centered_gaussian_window_cut(x::Array, min_x::Float64, max_x::Float64, n_σ::Real, center::Float64=0.0, n_bins_cut::Int=500, relative_cut::Float64=0.2, left::Bool=false)

    Cut out a single peak from the array x between min_x and max_x by fitting a truncated one-sided Gaussian and extrapolating a window cut with n_σ standard deviations. The center and side of the fit can be specified with left and center variable.

    Returns

    * `low_cut`: lower edge of the cut peak
    * `high_cut`: upper edge of the cut peak
    * `center`: center of the peak
    * `σ`: standard deviation of the Gaussian
    * `low_cut_fit`: lower edge of the cut peak from the fit
    * `high_cut_fit`: upper edge of the cut peak from the fit
    * `err`: error of the fit parameters
    source
    LegendSpecFits.get_continuum_surrival_fractionMethod
    get_continuum_surrival_fraction(aoe:::Vector{<:Unitful.RealOrRealQuantity}, e::Vector{<:T}, center::T, window::T, aoe_cut::Unitful.RealOrRealQuantity,; sigma_high_sided::Unitful.RealOrRealQuantity=NaN) where T<:Unitful.Energy{<:Real}

    Get the surrival fraction of a continuum after a AoE cut value aoe_cut for a given center and window size.

    Returns

    • center: Center of the continuum
    • window: Window size
    • n_before: Number of counts before the cut
    • n_after: Number of counts after the cut
    • sf: Surrival fraction
    source
    LegendSpecFits.get_mc_value_shapesMethod
    get_mc_value_shapes(v::NamedTuple, v_err::Matrix, n::Integer)

    Generate n random samples of fit parameters using their respective best-fit values v and covariance matrix v_err

    source
    LegendSpecFits.get_mc_value_shapesMethod
    get_mc_value_shapes(v::NamedTuple, v_err::NamedTuple, n::Integer)

    Return a NamedTuple with the same fields as v and v_err but with Normal distributions for each field.

    source
    LegendSpecFits.get_n_after_aoe_cutMethod
    get_n_after_aoe_cut(aoe_cut::T, aoe::Array{T}, e::Array{T}, peak::T, window::Array{T}, bin_width::T, result_before::NamedTuple, peakstats::NamedTuple; uncertainty=true) where T<:Real

    Get the number of counts after a AoE cut value aoe_cut for a given peak and window size while performing a peak fit with fixed position. The number of counts is determined by fitting the peak with a pseudo prior for the peak position.

    Returns

    • n: Number of counts after the cut
    source
    LegendSpecFits.get_peak_fwhm_th228Function
    get_peak_fwhm_th228(v_ml::NamedTuple, v_ml_err::NamedTuple)

    Get the FWHM of a peak from the fit parameters while performing a MC error propagation.

    Returns

    * `fwhm`: the FWHM of the peak
    * `fwhm_err`: the uncertainty of the FWHM of the peak
    source
    LegendSpecFits.get_peak_surrival_fractionMethod
    get_peak_surrival_fraction(aoe::Array{T}, e::Array{T}, peak::T, window::Array{T}, aoe_cut::T,; uncertainty::Bool=true, low_e_tail::Bool=true) where T<:Real

    Get the surrival fraction of a peak after a AoE cut value aoe_cut for a given peak and window size whiile performing a peak fit with fixed position.

    Returns

    • peak: Peak position
    • n_before: Number of counts before the cut
    • n_after: Number of counts after the cut
    • sf: Surrival fraction
    • err: Uncertainties
    source
    LegendSpecFits.get_peakhists_th228Method
    get_peakhists_th228(e::Array, th228_lines::Array, window_sizes::Array, e_unit::String="keV", proxy_binning_peak::Float64=2103.5, proxy_binning_peak_window::Float64=10.0)

    Create histograms around the calibration lines and return the histograms and the peak statistics.

    Returns

    * `peakhists`: array of histograms around the calibration lines
    * `peakstats`: array of statistics for the calibration line fits
    source
    LegendSpecFits.get_peaks_surrival_fractionsMethod
    get_peaks_surrival_fractions(aoe::Array{T}, e::Array{T}, peaks::Array{T}, peak_names::Array{Symbol}, windows::Array{T}, aoe_cut::T,; uncertainty=true) where T<:Real

    Get the surrival fraction of a peak after a AoE cut value aoe_cut for a given peak and window size while performing a peak fit with fixed position.

    Return

    • result: Dict of results for each peak
    • report: Dict of reports for each peak
    source
    LegendSpecFits.get_residualsMethod
    residuals(f_fit::Base.Callable, h::Histogram{<:Real,1},v_ml::NamedTuple)

    Calculate bin-wise residuals and normalized residuals. Calcualte bin-wise p-value based on poisson distribution for each bin.

    Input:

    • f_fitfunction handle of fit function (peakshape)
    • h histogram of data
    • v_ml best-fit parameters

    Returns:

    • residuals difference: model - data (histogram bin count)
    • residuals_norm normalized residuals: model - data / sqrt(model)
    • p_value_binwise p-value for each bin based on poisson distribution
    source
    LegendSpecFits.get_sf_after_aoe_cutMethod
    get_sf_after_aoe_cut(aoe_cut::T, aoe::Array{T}, e::Array{T}, peak::T, window::Array{T}, bin_width::T, result_before::NamedTuple; uncertainty=true) where T<:Real

    Get the survival fraction after a AoE cut value aoe_cut for a given peak and window size from a combined fit to the survived and cut histograms.

    Returns

    • sf: Survival fraction after the cut
    source
    LegendSpecFits.hist_loglikeMethod
    hist_loglike(f_fit::Base.Callable, h::Histogram{<:Real,1})

    Calculate the Poisson log-likelihood of a fit function f_fit(x) and a histogram h. f_fit must accept all values x on the horizontal axis of the histogram.

    Currently uses a simple midpoint-rule integration of f_fit over the bins of h.

    source
    LegendSpecFits.lowEtail_peakshapeMethod
    lowEtail_peakshape(
        x::Real, μ::Real, σ::Real, n::Real,
        skew_fraction::Real, skew_width::Real,
    )

    Describes the low-E signal tail part of the shape of a typical gamma peak in a detector.

    source
    LegendSpecFits.nearestSPDMethod
    nearestSPD(A::Matrix{<:Real})

    Returns the nearest positive definite matrix to A Calculation is based on matrix factorization techniques described in https://www.sciencedirect.com/science/article/pii/0024379588902236

    source
    LegendSpecFits.p_valueMethod
    p_value(f_fit::Base.Callable, h::Histogram{<:Real,1},v_ml::NamedTuple)

    calculate p-value based on least-squares, assuming poisson uncertainty baseline method to get goodness-of-fit (gof)

    input:

    • f_fitfunction handle of fit function (peakshape)
    • h histogram of data
    • v_ml best-fit parameters

    returns:

    • pval p-value of chi2 test
    • chi2 chi2 value
    • dof degrees of freedom
    source
    LegendSpecFits.p_value_MCMethod
    p_value_MC(f_fit::Base.Callable, h::Histogram{<:Real,1},ps::NamedTuple{(:peak_pos, :peak_fwhm, :peak_sigma, :peak_counts, :mean_background)},v_ml::NamedTuple,;n_samples::Int64=1000)

    alternative p-value calculation via Monte Carlo sampling. Warning: computational more expensive than pvaule() and pvalue_LogLikeRatio()

    Input:

    • f_fitfunction handle of fit function (peakshape)
    • h histogram of data
    • ps best-fit parameters
    • v_ml best-fit parameters
    • n_samples number of samples

    Performed Steps:

    • Create n_samples randomized histograms. For each bin, samples are drawn from a Poisson distribution with λ = model peak shape (best-fit parameter)
    • Each sample histogram is fit using the model function f_fit
    • For each sample fit, the max. loglikelihood fit is calculated

    Returns

    • % p value –> comparison of sample max. loglikelihood and max. loglikelihood of best-fit
    source
    LegendSpecFits.prepare_dep_peakhistMethod
    prepare_dep_peakhist(e::Array{T}, dep::T,; relative_cut::T=0.5, n_bins_cut::Int=500) where T<:Real

    Prepare an array of uncalibrated DEP energies for parameter extraction and calibration.

    Returns

    • result: Result of the initial fit
    • report: Report of the initial fit
    source
    LegendSpecFits.pulser_cal_qcMethod
    pulser_cal_qc(data, pulser_config; n_pulser_identified=100)

    Perform simple QC cuts on the data and return the data for energy calibration.

    Returns

    - pulser_idx: indices of the pulser events
    source
    LegendSpecFits.signal_peakshapeMethod
    signal_peakshape(
        x::Real, μ::Real, σ::Real, n::Real,
        skew_fraction::Real, skew_width::Real,
    )

    Describes the signal part of the shape of a typical gamma peak in a detector.

    source
    LegendSpecFits.simple_calibrationFunction
    simple_calibration(e_uncal::Array, th228_lines::Array, window_size::Float64=25.0, n_bins::Int=15000, calib_type::String="th228")

    Perform a simple calibration for the uncalibrated energy array e_uncal using the calibration type calib_type and the calibration lines th228_lines. The window size is the size of the window around the calibration line to use for the calibration. The number of bins is the number of bins to use for the histogram.

    Returns * h_calsimple: histogram of the calibrated energy array * h_uncal: histogram of the uncalibrated energy array * c: calibration factor * fep_guess: estimated full energy peak (FEP) * peakhists: array of histograms around the calibration lines * peakstats: array of statistics for the calibration line fits

    source
    LegendSpecFits.step_gaussMethod
    step_gauss(x::Real, μ::Real, σ::Real)

    Evaluates the convulution of a Heaviside step function and the PDF of Normal(μ, σ) at x.

    The result does not correspond to a PDF as it is not normalizable.

    source
    LegendSpecFits.tuple_to_arrayMethod
    tuple_to_array(nt::NamedTuple, fields::Vector{Symbol})

    Return an array with the values of the fields in nt in the order given by fields.

    source
    LegendSpecFits.weibull_from_mxFunction
    weibull_from_mx(m::Real, x::Real, p_x::Real = 0.6827)::Weibull

    Construct a Weibull distribution with a given median m and a given p_x-quantile x.

    Useful to construct priors for positive quantities.

    source