Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Constructor using Dict{Symbol, Key} from UnROOT's hist #32

Closed
sgnoohc opened this issue Aug 26, 2021 · 5 comments
Closed

Constructor using Dict{Symbol, Key} from UnROOT's hist #32

sgnoohc opened this issue Aug 26, 2021 · 5 comments

Comments

@sgnoohc
Copy link

sgnoohc commented Aug 26, 2021

Below is an example of how one might turn a UnROOT's hist object into a Hist1D.

function make_fhist1d(th1::Dict{Symbol, Any})
    xmin = th1[:fXaxis_fXmin]
    xmax = th1[:fXaxis_fXmax]
    nbins = th1[:fXaxis_fNbins]
    bins = range(xmin, xmax, length=nbins+1)
    h = Hist1D(Float64; bins=bins)
    h.hist.weights .= th1[:fN][2:end-1]
    h.sumw2 .= th1[:fSumw2][2:end-1]
    return h
end

I am not sure what would be the best way or best place to incorporate stuff like this....

@sgnoohc
Copy link
Author

sgnoohc commented Aug 26, 2021

I don't know how to get the bin edges directly rather than computing from xmin xmax and nbins.
I couldn't find the key in th1 from UnROOT.
(but must be there somewhere..?)

@Moelf
Copy link
Owner

Moelf commented Aug 26, 2021

I think we could use a method that returns ( edges counts and sumw2 ) in UnROOT.jl, since it's easy to construct histogram from those in either StatsBase or FHist and it doesn't bloat the dependency.

the individuals bin edges are indeed not stored when they are uniform, your observation is correct I think

@aminnj
Copy link
Collaborator

aminnj commented Aug 26, 2021

I agree with @Moelf. About the irregular binning, this is how to get it:

julia> using UnROOT, FHist

julia> f = ROOTFile("/Users/namin/.julia/dev/UnROOT/test/samples/histograms1d2d.root")

julia> th = f["myTH1D_nonuniform"];

julia> th[:fXaxis_fXmin], th[:fXaxis_fXmax], th[:fXaxis_fNbins], th[:fXaxis_fXbins]
(-2.0, 2.0, 2, [-2.0, 1.0, 2.0])

julia> th = f["myTH1D"];

julia> th[:fXaxis_fXmin], th[:fXaxis_fXmax], th[:fXaxis_fNbins], th[:fXaxis_fXbins]
(-2.0, 2.0, 2, Any[])

That is, :fXaxis_fXbins is filled if the histogram has irregular binning, otherwise min/max/nbins is sufficient to construct the edges with range(min,max,length=nbins+1)

@aminnj
Copy link
Collaborator

aminnj commented Aug 26, 2021

Another thing to be careful about if we want a method that returns edges, counts, sumw2 in UnROOT.jl: what do we do with overflow bins? If we return the full arrays as-is, then the user still has to go through and do counts[2:end-1] and sumw2[2:end-1]

@aminnj
Copy link
Collaborator

aminnj commented Aug 26, 2021

based on @sgnoohc's example, but only half-tested:

function make_fhist1d(th1::Dict{Symbol, Any}; overflow=false)
    xmin = th1[:fXaxis_fXmin]
    xmax = th1[:fXaxis_fXmax]
    nbins = th1[:fXaxis_fNbins]
    bins = isempty(th1[:fXaxis_fXbins]) ? range(xmin, xmax, length=nbins+1) : th1[:fXaxis_fXbins];
    h = Hist1D(Float64; bins=bins, overflow=overflow) # so if we push more later, overflow is respected
    counts = th1[:fN]
    sumw2 = th1[:fSumw2]
    if overflow
        counts[2] += counts[1]
        counts[end-1] += counts[end]
        sumw2[2] += sumw2[1]
        sumw2[end-1] += sumw2[end]
    end
    h.hist.weights .= counts[2:end-1]
    h.sumw2 .= sumw2[2:end-1]
    return h
end

and then 2D is to be done...

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

No branches or pull requests

3 participants