FHist.ATLASTHEME
— ConstantExample
with_theme(ATLASTHEME) do
h1 = Hist1D(randn(10^4))
hist(h1; label="atlas style histogram")
end
FHist.BinEdges
— TypeBinEdges <: AbstractVector{Float64}
This type implements a vector-like data structure to be used for histogram bin edges, it can handle both uniform and non-uniform binnings in a single type to reduce the amount of parametric types. It would switch to O(1) `searchsortedlast` if the binning is uniform.
Due to the usage of Float64, bin edges shouldn't contain element with absolute value larger than 9007199254740992, which is the maxintfloat(Float64)
.
FHist.Hist1D
— MethodTo make an empty histogram
use the all-keyword-arguments constructor:
Hist1D(;
counttype=Float64,
binedges::E
bincounts = zeros(counttype, length.(_to_tuple(binedges)) .- 1),
sumw2 = zero(bincounts),
nentries = 0
overflow::Bool = false
) where {E<:NTuple{1,Any}}
Everything other than binedges
are optional (infered from binedges
).
To make an histogram given data (and weights etc.)
use the a positional argument for data and keyword-arguments for the rest:
Hist1D(array::E;
counttype=Float64,
binedges=nothing,
weights=nothing,
nbins=nothing,
overflow=false)
) where {E<:NTuple{1,Any}}
Everything other than data (array
) is optional (infered from data).
FHist.Hist2D
— MethodTo make an empty histogram
use the all-keyword-arguments constructor:
Hist2D(;
counttype=Float64,
binedges::E
bincounts = zeros(counttype, length.(_to_tuple(binedges)) .- 1),
sumw2 = zero(bincounts),
nentries = 0
overflow::Bool = false
) where {E<:NTuple{2,Any}}
Everything other than binedges
are optional (infered from binedges
).
To make an histogram given data (and weights etc.)
use the a positional argument for data and keyword-arguments for the rest:
Hist2D(array::E;
counttype=Float64,
binedges=nothing,
weights=nothing,
nbins=nothing,
overflow=false)
) where {E<:NTuple{2,Any}}
Everything other than data (array
) is optional (infered from data).
FHist.Hist3D
— MethodTo make an empty histogram
use the all-keyword-arguments constructor:
Hist3D(;
counttype=Float64,
binedges::E
bincounts = zeros(counttype, length.(_to_tuple(binedges)) .- 1),
sumw2 = zero(bincounts),
nentries = 0
overflow::Bool = false
) where {E<:NTuple{3,Any}}
Everything other than binedges
are optional (infered from binedges
).
To make an histogram given data (and weights etc.)
use the a positional argument for data and keyword-arguments for the rest:
Hist3D(array::E;
counttype=Float64,
binedges=nothing,
weights=nothing,
nbins=nothing,
overflow=false)
) where {E<:NTuple{3,Any}}
Everything other than data (array
) is optional (infered from data).
FHist._factor
— Methodfunction _factor(n::Integer)
Helper function to calculate the prime factors of a given integer.
FHist.atomic_push!
— Functionpush!(h::Hist2D, valx::Real, valy::Real, wgt::Real=1)
atomic_push!(h::Hist2D, valx::Real, valy::Real, wgt::Real=1)
Adding one value at a time into histogram. sumw2
(sum of weights^2) accumulates wgt^2
with a default weight of 1. atomic_push!
is a slower version of push!
that is thread-safe.
FHist.atomic_push!
— Functionpush!(h::Hist3D, valx::Real, valy::Real, wgt::Real=1)
atomic_push!(h::Hist3D, valx::Real, valy::Real, wgt::Real=1)
Adding one value at a time into histogram. sumw2
(sum of weights^2) accumulates wgt^2
with a default weight of 1. atomic_push!
is a slower version of push!
that is thread-safe.
FHist.atomic_push!
— Functionpush!(h::Hist1D, val::Real, wgt::Real=1)
atomic_push!(h::Hist1D, val::Real, wgt::Real=1)
Adding one value at a time into histogram. sumw2
(sum of weights^2) accumulates wgt^2
with a default weight of 1. atomic_push!
is a slower version of push!
that is thread-safe.
N.B. To append multiple values at once, use broadcasting via push!.(h, [-3.0, -2.9, -2.8])
or push!.(h, [-3.0, -2.9, -2.8], 2.0)
FHist.bayes_rebin_edges
— Methodbayes_rebin_edges(h::Hist1D; prior=BayesHistogram.Geometric(0.995))
Find optimal bin edges for a histogram using Bayesian rebinning algorithm. This function only find edges, it doesn't return a new histogram.
For possible priors, see BayesHistogram.jl
.
FHist.bincenters
— Methodbincenters(h::Hist1D)
Get the bin centers of the histogram
For 1D histogram, it returns just a vector. For others, it returns a tuple of vectors.
FHist.bincenters
— Methodbincenters(h::Hist2D)
Get the bin centers of the histogram
For 1D histogram, it returns just a vector. For others, it returns a tuple of vectors.
FHist.bincenters
— Methodbincenters(h::Hist3D)
Get the bin centers of the histogram
For 1D histogram, it returns just a vector. For others, it returns a tuple of vectors.
FHist.bincounts
— Method bincounts(h::Hist1D)
Get the bin counts (weights) of the histogram.
FHist.bincounts
— Method bincounts(h::Hist2D)
Get the bin counts (weights) of the histogram.
FHist.bincounts
— Method bincounts(h::Hist3D)
Get the bin counts (weights) of the histogram.
FHist.binedges
— Methodbinedges(h)
Get the bin edges of the histogram
For 1D histogram, it returns just a vector. For others, it returns a tuple of vectors. If you need a tuple of vectors, use h.binedges
at your own risk.
FHist.binedges
— Methodbinedges(h)
Get the bin edges of the histogram
For 1D histogram, it returns just a vector. For others, it returns a tuple of vectors. If you need a tuple of vectors, use h.binedges
at your own risk.
FHist.binedges
— Methodbinedges(h)
Get the bin edges of the histogram
For 1D histogram, it returns just a vector. For others, it returns a tuple of vectors. If you need a tuple of vectors, use h.binedges
at your own risk.
FHist.binerrors
— Methodbinerrors(f=sqrt, h)
Get the error (uncertainty) of each bin. By default, calls sqrt
on sumw2(h)
bin by bin as an approximation.
FHist.binerrors
— Methodbinerrors(f=sqrt, h)
Get the error (uncertainty) of each bin. By default, calls sqrt
on sumw2(h)
bin by bin as an approximation.
FHist.binerrors
— Methodbinerrors(f=sqrt, h)
Get the error (uncertainty) of each bin. By default, calls sqrt
on sumw2(h)
bin by bin as an approximation.
FHist.cumulative
— Methodcumulative(h::Hist1D; forward=true)
Create a cumulative histogram. If forward
, start summing from the left.
FHist.effective_entries
— Methodeffective_entries(h) -> scalar
Get the number of effective entries for the entire histogram:
\[n_{eff} = \frac{(\sum Weights )^2}{(\sum Weight^2 )}\]
This is also equivalent to integral(hist)^2 / sum(sumw2(hist))
, this is the same as TH1::GetEffectiveEntries()
FHist.effective_entries
— Methodeffective_entries(h) -> scalar
Get the number of effective entries for the entire histogram:
\[n_{eff} = \frac{(\sum Weights )^2}{(\sum Weight^2 )}\]
This is also equivalent to integral(hist)^2 / sum(sumw2(hist))
, this is the same as TH1::GetEffectiveEntries()
FHist.effective_entries
— Methodeffective_entries(h) -> scalar
Get the number of effective entries for the entire histogram:
\[n_{eff} = \frac{(\sum Weights )^2}{(\sum Weight^2 )}\]
This is also equivalent to integral(hist)^2 / sum(sumw2(hist))
, this is the same as TH1::GetEffectiveEntries()
FHist.hists_to_bars
— Methodhists_to_bars(hist1ds)
Given a vector of Hist1D, return edges
(xs), heights
(ys), and grps
(for grouping) that is useful for plotting stacked histogram.
FHist.integral
— Methodintegral(h; width=false)
Get the integral a histogram; width
means multiply each bincount by their bin width when calculating the integral.
width
keyword argument only works with 1D histogram at the moment.
Be aware of the approximation you make when using width=true
with histogram with overflow bins, the overflow bins (i.e. the left/right most bins) width will be taken "as is".
FHist.lookup
— Methodlookup(h::Hist1D, x)
For given x-axis value x
, find the corresponding bin and return the bin content. If a value is out of the histogram range, return missing
.
FHist.lookup
— Methodfunction lookup(h::Hist2D, x, y)
For given x-axis and y-axis value x
, y
, find the corresponding bin and return the bin content. If a value is out of the histogram range, return missing
.
FHist.lookup
— Methodfunction lookup(h::Hist3D, x, y, z)
For given x/y/z-axis value x
, y
, z
, find the corresponding bin and return the bin content. If a value is out of the histogram range, return missing
.
FHist.nbins
— Methodnbins(h::Hist1D)
Get the number of bins of a histogram.
FHist.nbins
— Methodnbins(h::Hist2D)
Get a 2-tuple of the number of x and y bins of a histogram.
FHist.nbins
— Methodnbins(h::Hist3D)
Get a 3-tuple of the number of x and y bins of a histogram.
FHist.nentries
— Methodnentries(h::Hist1D)
Get the number of times a histogram is filled (push!
ed)
FHist.nentries
— Methodnentries(h::Hist2D)
Get the number of times a histogram is filled (push!
ed)
FHist.nentries
— Methodnentries(h::Hist3D)
Get the number of times a histogram is filled (push!
ed)
FHist.profile
— Functionprofile(h::Hist2D, axis::Symbol=:x)
profile(axis::Symbol=:x) = h::Hist2D -> profile(h, axis)
Returns the axis
-profile of the 2D histogram by calculating the weighted mean over the other axis. profile(h, :x)
will return a Hist1D
with the y-axis edges of h
.
FHist.project
— Functionproject(h::Hist2D, axis::Symbol=:x)
project(axis::Symbol=:x) = h::Hist2D -> project(h, axis)
Computes the :x
(:y
) axis projection of the 2D histogram by summing over the y (x) axis. Returns a Hist1D
.
FHist.project
— Functionproject(h::Hist3D, axis::Symbol=:x)
project(axis::Symbol=:x) = h::Hist3D -> project(h, axis)
Computes the :x
/:y
/:z
axis projection of the 3D histogram by summing over the specified axis. Returns a Hist2D
.
FHist.rebin
— Functionrebin(h::Hist2D, nx::Int=1, ny::Int=nx)
rebin(nx::Int, ny::Int) = h::Hist2D -> rebin(h, nx, ny)
Merges nx
(ny
) consecutive bins into one along the x (y) axis by summing.
FHist.rebin
— Functionrebin(h::Hist1D, n::Int=1)
rebin(n::Int) = h::Hist1D -> rebin(h, n)
Merges n
consecutive bins into one. The returned histogram will have nbins(h)/n
bins.
FHist.restrict
— Functionrestrict(h::Hist2D, xlow=-Inf, xhigh=Inf, ylow=-Inf, yhigh=Inf)
restrict(xlow=-Inf, xhigh=Inf, ylow=-Inf, yhigh=Inf) = h::Hist2D -> restrict(h, xlow, xhigh, ylow, yhigh)
Returns a new histogram with a restricted x-axis. restrict(h, 0, 3)
(or h |> restrict(0, 3)
) will return a slice of h
where the bin centers are in [0, 3]
(inclusive).
FHist.restrict
— Functionrestrict(h::Hist1D, low=-Inf, high=Inf)
restrict(low=-Inf, high=Inf) = h::Hist1D -> restrict(h, low, high)
Returns a new histogram with a restricted x-axis. restrict(h, 0, 3)
(or h |> restrict(0, 3)
) will return a slice of h
where the bin centers are in [0, 3]
(inclusive).
FHist.sample
— Methodsample(h::Hist1D, n::Int=1)
Sample a histogram's with weights equal to bin count, n
times. The sampled values are the bins' lower edges.
FHist.significance
— Methodsignificance(signal, bkg) -> `(significance, error_on_significance)`
Calculate the significance of signal vs. bkg histograms, this function uses a more accurate algorithm than the naive S / sqrt(B)
Ref: https://cds.cern.ch/record/2736148/files/ATL-PHYS-PUB-2020-025.pdf
Example:
h1 = Hist1D(rand(1000); binedges = [0, 0.5])
h2 = Hist1D(rand(10000); binedges = [0, 0.5]);
julia> s1 = significance(h1,h2)
(6.738690967342175, 0.3042424717261312)
FHist.sumw2
— Methodsumw2(h)
Get the sum of weights squared of the histogram, it has the same shape as bincounts(h)
.
FHist.sumw2
— Methodsumw2(h)
Get the sum of weights squared of the histogram, it has the same shape as bincounts(h)
.
FHist.sumw2
— Methodsumw2(h)
Get the sum of weights squared of the histogram, it has the same shape as bincounts(h)
.
FHist.transpose
— Methodtranspose(h::Hist2D)
Reverses the x and y axes.
FHist.valid_rebin_values
— Methodvalid_rebin_values(h::Union{Hist1D, Hist2D, Hist3D})
Calculates the legal values for rebinning, essentially the prime factors of the number of bins. For a 1D histogram, a Set
of numbers is return, for higher dimensional histograms a Vector{Set}
for each dimension.
LinearAlgebra.normalize
— Methodnormalize(h::Hist1D; width=true)
Create a normalized histogram via division by integral(h)
, when width==true
, the resultant histogram has area under the curve equals 1.
Implicit approximation is made when using width=true
with histograms that have overflow bins: the overflow data lives inthe left/right most bins and the bin width is taken "as is".
LinearAlgebra.normalize
— Methodnormalize(h::Hist2D)
Create a normalized histogram via division by integral(h)
.
LinearAlgebra.normalize
— Methodnormalize(h::Hist3D)
Create a normalized histogram via division by integral(h)
.
Statistics.mean
— MethodStatistics.mean(h)
Statistics.std(h)
Statistics.median(h)
Statistics.quantile(h::Hist1D, p)
Compute statistical quantities based on the bin centers weighted by the bin counts.
When the histogram is Hist2D
, return tuple instead, e.g (mean(project(h, :x)), mean(project(h, :y)))
etc.