LegendSpecFits
Also see LegendSpecFits on GitHub and the full package documentation.
Modules
Types and constants
Functions and macros
LegendSpecFits._get_model_counts
LegendSpecFits._prepare_data
LegendSpecFits.aoe_compton_background_peakshape
LegendSpecFits.aoe_compton_peakshape
LegendSpecFits.aoe_compton_signal_peakshape
LegendSpecFits.array_to_tuple
LegendSpecFits.autocal_energy
LegendSpecFits.background_peakshape
LegendSpecFits.baseline_qc
LegendSpecFits.chi2fit
LegendSpecFits.ctc_energy
LegendSpecFits.cut_single_peak
LegendSpecFits.estimate_combined_peak_stats
LegendSpecFits.estimate_fwhm
LegendSpecFits.estimate_single_peak_stats
LegendSpecFits.estimate_single_peak_stats_psd
LegendSpecFits.ex_gauss_pdf
LegendSpecFits.ex_step_gauss
LegendSpecFits.expand_vars
LegendSpecFits.exponential_decay
LegendSpecFits.f_optimize_ctc
LegendSpecFits.fit_aoe_compton
LegendSpecFits.fit_aoe_corrections
LegendSpecFits.fit_calibration
LegendSpecFits.fit_enc_sigmas
LegendSpecFits.fit_fwhm
LegendSpecFits.fit_fwhm_ft_fep
LegendSpecFits.fit_half_centered_trunc_gauss
LegendSpecFits.fit_half_trunc_gauss
LegendSpecFits.fit_peaks
LegendSpecFits.fit_peaks_combined
LegendSpecFits.fit_sg_wl
LegendSpecFits.fit_single_aoe_compton
LegendSpecFits.fit_single_peak_th228
LegendSpecFits.fit_single_trunc_gauss
LegendSpecFits.fit_subpeaks_th228
LegendSpecFits.gamma_peakshape
LegendSpecFits.gauss_pdf
LegendSpecFits.generate_aoe_compton_bands
LegendSpecFits.get_aoe_cut
LegendSpecFits.get_centered_gaussian_window_cut
LegendSpecFits.get_continuum_surrival_fraction
LegendSpecFits.get_distribution_transform
LegendSpecFits.get_friedman_diaconis_bin_width
LegendSpecFits.get_mc_value_shapes
LegendSpecFits.get_mc_value_shapes
LegendSpecFits.get_n_after_aoe_cut
LegendSpecFits.get_number_of_bins
LegendSpecFits.get_peak_fwhm_th228
LegendSpecFits.get_peak_surrival_fraction
LegendSpecFits.get_peakhists_th228
LegendSpecFits.get_peaks_surrival_fractions
LegendSpecFits.get_residuals
LegendSpecFits.get_sf_after_aoe_cut
LegendSpecFits.hist_loglike
LegendSpecFits.linear_function
LegendSpecFits.lowEtail_peakshape
LegendSpecFits.nearestSPD
LegendSpecFits.p_value
LegendSpecFits.p_value_LogLikeRatio
LegendSpecFits.p_value_MC
LegendSpecFits.prepare_dep_peakhist
LegendSpecFits.pulser_cal_qc
LegendSpecFits.signal_peakshape
LegendSpecFits.simple_calibration
LegendSpecFits.step_gauss
LegendSpecFits.subhist
LegendSpecFits.tuple_to_array
LegendSpecFits.weibull_from_mx
Documentation
LegendSpecFits.LegendSpecFits
— ModuleLegendSpecFits
Template for Julia packages.
LegendSpecFits._get_model_counts
— Method_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
LegendSpecFits._prepare_data
— Method_prepare_data(h::Histogram{<:Real,1})
aux. function to convert histogram data into bin edges, bin width and bin counts
LegendSpecFits.aoe_compton_background_peakshape
— Methodaoe_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.
LegendSpecFits.aoe_compton_peakshape
— Methodaoe_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.
LegendSpecFits.aoe_compton_signal_peakshape
— Methodaoe_compton_signal_peakshape(
x::Real, μ::Real, σ::Real, n::Real
)
Describes the signal shape of a typical A/E Compton peak in a detector as a gaussian SSE peak.
LegendSpecFits.array_to_tuple
— Methodarray_to_tuple(a::AbstractArray, as_nt::NamedTuple)
Return a NamedTuple
with the values of a
in the order given by fieldnames(as_nt)
.
LegendSpecFits.autocal_energy
— Methodautocal_energy(E_raw::AbstractArray{<:Real})
Compute an energy calibration from raw reconstructed energy deposition values.
LegendSpecFits.background_peakshape
— Methodbackground_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.
LegendSpecFits.baseline_qc
— Methodbaseline_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 variablereport
: Namedtuple containing plotable objects.
LegendSpecFits.chi2fit
— Methodfit_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:
LegendSpecFits.ctc_energy
— Functionctc_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
LegendSpecFits.cut_single_peak
— Methodcut_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
LegendSpecFits.estimate_combined_peak_stats
— Methodestimate_combined_peak_stats(peakstats::StructArray,; calib_type::Symbol=:th228)
Estimate the peak position, FWHM, sigma, counts and background of a peak from a histogram.
LegendSpecFits.estimate_fwhm
— Methodestimate_fwhm(v::NamedTuple, v_err::NamedTuple)
Get the FWHM of a peak from the fit parameters.
Returns
* `fwhm`: the FWHM of the peak
LegendSpecFits.estimate_single_peak_stats
— Methodestimate_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
LegendSpecFits.estimate_single_peak_stats_psd
— Methodestimate_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
LegendSpecFits.ex_gauss_pdf
— Methodex_gauss_pdf(x::Real, μ::Real, σ::Real, θ::Real)
The PDF of an Exponentially modified Gaussian distribution with Gaussian parameters μ
, σ
and exponential scale θ
at x
.
It is the PDF of the distribution that descibes the random process rand(Normal(μ, σ)) + rand(Exponential(θ))
.
LegendSpecFits.ex_step_gauss
— Methodex_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
.
LegendSpecFits.expand_vars
— Methodexpand_vars(v::NamedTuple)::StructArray
Expand all fields in v
(scalars or arrays) to same array size and return a StructArray
.
LegendSpecFits.exponential_decay
— Methodexponential_decay(x::Real, amplitude::Real, decay::Real, offset::Real)
Evaluates an exponential decay function at x
with parameters amplitude
, decay
and offset
.
LegendSpecFits.f_optimize_ctc
— Methodf_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
LegendSpecFits.fit_aoe_compton
— Methodfit_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
LegendSpecFits.fit_aoe_corrections
— Methodfitaoecorrections(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 valuesf_μ_scs
: Fit function for the mean valuesσ_scs
: Fit result for the sigma valuesf_σ_scs
: Fit function for the sigma values
LegendSpecFits.fit_calibration
— Methodfit_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`:
LegendSpecFits.fit_enc_sigmas
— Methodfit_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 inenc_grid_rt
enc_grid_rt
: 1D array of RT values for which the ENC values inenc_grid
are calculatedmin_enc
: minimum ENC value to consider for the fitmax_enc
: maximum ENC value to consider for the fitnbins
: number of bins to use for the histogram of ENC valuesrel_cut_fit
: relative cut value to use for the fit
Returns
rt
: optimal RT valuemin_enc
: corresponding ENC value
LegendSpecFits.fit_fwhm
— MethodfitFWHM(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
LegendSpecFits.fit_fwhm_ft_fep
— Methodfitfwhmftfep(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 ine_grid_ft
e_grid_ft
: 1D array of FT values for which the FWHM values ine_grid
are calculatedrt
: RT value for which the FWHM values ine_grid
are calculatedmin_e
: minimum energy value to consider for the fitmax_e
: maximum energy value to consider for the fitnbins
: number of bins to use for the histogram of energy valuesrel_cut_fit
: relative cut value to use for the fit
Returns
ft
: optimal FT valuemin_fwhm
: corresponding FWHM value
LegendSpecFits.fit_half_centered_trunc_gauss
— Methodfit_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
LegendSpecFits.fit_half_trunc_gauss
— Methodfit_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
LegendSpecFits.fit_peaks
— Methodfit_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
LegendSpecFits.fit_peaks_combined
— Methodfit_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.
LegendSpecFits.fit_sg_wl
— Methodfit_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 dataa_grid_wl_sg
: range of window lengths to sweep throughoptimization_config
: configuration dictionary
Returns
result
: optimal window length and corresponding survival fractionreport
: report with all window lengths and survival fractions
LegendSpecFits.fit_single_aoe_compton
— Methodfit_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
LegendSpecFits.fit_single_peak_th228
— Methodfit_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
LegendSpecFits.fit_single_trunc_gauss
— Methodfit_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 result
with: *
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
LegendSpecFits.fit_subpeaks_th228
— Methodfit_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
LegendSpecFits.gamma_peakshape
— Methodgamma_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.
LegendSpecFits.gauss_pdf
— MethodLegendSpecFits.gauss_pdf(x::Real, μ::Real, σ::Real)
Equivalent to pdf(Normal(μ, σ), x)
LegendSpecFits.generate_aoe_compton_bands
— Methodgenerate_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
LegendSpecFits.get_aoe_cut
— Methodget_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 valuen0
: Number of counts before the cutnsf
: Number of counts after the cut
LegendSpecFits.get_centered_gaussian_window_cut
— Methodget_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
LegendSpecFits.get_continuum_surrival_fraction
— Methodget_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 continuumwindow
: Window sizen_before
: Number of counts before the cutn_after
: Number of counts after the cutsf
: Surrival fraction
LegendSpecFits.get_distribution_transform
— Functionget_distribution_transform(d::Distribution, pprior::Prior)
Return a DistributionTransform
for the given Distribution
and Prior
.
LegendSpecFits.get_friedman_diaconis_bin_width
— Functionget_friedman_diaconis_bin_width(x::AbstractArray)
Return the bin width for the given data x
using the Friedman-Diaconis rule.
LegendSpecFits.get_mc_value_shapes
— Methodget_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
LegendSpecFits.get_mc_value_shapes
— Methodget_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.
LegendSpecFits.get_n_after_aoe_cut
— Methodget_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
LegendSpecFits.get_number_of_bins
— Methodget_number_of_bins(x::AbstractArray,; method::Symbol=:sqrt)
Return the number of bins for the given data x
using the given method.
LegendSpecFits.get_peak_fwhm_th228
— Functionget_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
LegendSpecFits.get_peak_surrival_fraction
— Methodget_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 positionn_before
: Number of counts before the cutn_after
: Number of counts after the cutsf
: Surrival fractionerr
: Uncertainties
LegendSpecFits.get_peakhists_th228
— Methodget_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
LegendSpecFits.get_peaks_surrival_fractions
— Methodget_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 peakreport
: Dict of reports for each peak
LegendSpecFits.get_residuals
— Methodresiduals(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_fit
function handle of fit function (peakshape)h
histogram of datav_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
LegendSpecFits.get_sf_after_aoe_cut
— Methodget_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
LegendSpecFits.hist_loglike
— Methodhist_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
.
LegendSpecFits.linear_function
— Methodlinear_function(x::Real, slope::Real, intercept::Real)
Evaluates a linear function at x
with parameters slope
and intercept
.
LegendSpecFits.lowEtail_peakshape
— MethodlowEtail_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.
LegendSpecFits.nearestSPD
— MethodnearestSPD(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
LegendSpecFits.p_value
— Methodp_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_fit
function handle of fit function (peakshape)h
histogram of datav_ml
best-fit parameters
returns:
pval
p-value of chi2 testchi2
chi2 valuedof
degrees of freedom
LegendSpecFits.p_value_LogLikeRatio
— Methodp_value_LogLikeRatio(f_fit::Base.Callable, h::Histogram{<:Real,1},v_ml::NamedTuple)
Alternative p-value via loglikelihood ratio
LegendSpecFits.p_value_MC
— Methodp_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_fit
function handle of fit function (peakshape)h
histogram of dataps
best-fit parametersv_ml
best-fit parametersn_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
LegendSpecFits.prepare_dep_peakhist
— Methodprepare_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 fitreport
: Report of the initial fit
LegendSpecFits.pulser_cal_qc
— Methodpulser_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
LegendSpecFits.signal_peakshape
— Methodsignal_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.
LegendSpecFits.simple_calibration
— Functionsimple_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
LegendSpecFits.step_gauss
— Methodstep_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.
LegendSpecFits.subhist
— Methodsubhist(h::Histogram, r::Tuple{<:Real,<:Real})
Return a new Histogram
with the bins in the range r
.
LegendSpecFits.tuple_to_array
— Methodtuple_to_array(nt::NamedTuple, fields::Vector{Symbol})
Return an array with the values of the fields in nt
in the order given by fields
.
LegendSpecFits.weibull_from_mx
— Functionweibull_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.