Skip to content

Commit

Permalink
Add TH3 support (#337)
Browse files Browse the repository at this point in the history
* Add TH3 support

* Remove debug @show

* Add TH3F sample file and script

* Improve TH3F sample

* Fix hist readout and add tests

* Allow returning missing for parents

* Add tests for TH3F

* More bootstrapping for ridiculous class inheritance madness

(aka "why do we need Latex stuff to read a 1D histo??!!??")
  • Loading branch information
tamasgal authored May 6, 2024
1 parent b94a2f8 commit 3aad70a
Show file tree
Hide file tree
Showing 5 changed files with 264 additions and 1 deletion.
194 changes: 194 additions & 0 deletions src/bootstrap.jl
Original file line number Diff line number Diff line change
Expand Up @@ -157,6 +157,9 @@ struct TH2_5 <: TH2 end
function readfields!(io, fields, T::Type{TH2_4}) end
function readfields!(io, fields, T::Type{TH2_5}) end

abstract type TH3 <: ROOTStreamedObject end
struct TH3_6 <: TH3 end
function readfields!(io, fields, T::Type{TH3_6}) end

@with_kw struct ROOT_3a3a_TIOFeatures <: ROOTStreamedObject
fIOBits
Expand Down Expand Up @@ -944,8 +947,10 @@ end

TH1F(io, tkey::TKey, refs) = TH(io, tkey, refs)
TH2F(io, tkey::TKey, refs) = TH(io, tkey, refs)
TH3F(io, tkey::TKey, refs) = TH(io, tkey, refs)
TH1D(io, tkey::TKey, refs) = TH(io, tkey, refs)
TH2D(io, tkey::TKey, refs) = TH(io, tkey, refs)
TH3D(io, tkey::TKey, refs) = TH(io, tkey, refs)

"""
TH(io, tkey::TKey, refs)
Expand All @@ -964,6 +969,11 @@ function TH(io, tkey::TKey, refs)
stream!(io, fields, TH2, check=false)
end

is3d = startswith(tkey.fClassName, "TH3")
if is3d
stream!(io, fields, TH3, check=false)
end

stream!(io, fields, TH1, check=false)
stream!(io, fields, TNamed)
stream!(io, fields, TAttLine)
Expand Down Expand Up @@ -1007,6 +1017,14 @@ function TH(io, tkey::TKey, refs)
end
end


if is3d
skip(io, 6) # TAtt3D, not used yet, consist of two fields: cnt (u4), vers (u2)
for symb in [:fTsumwy, :fTsumwy2, :fTsumwxy, :fTsumwz, :fTsumwz2, :fTsumwxz, :fTsumwyz]
fields[symb] = readtype(io, Float64)
end
end

arraytype = endswith(tkey.fClassName, 'F') ? TArrayF : TArrayD
fields[:fN] = readtype(io, arraytype)
fields
Expand Down Expand Up @@ -1221,3 +1239,179 @@ function readfields!(c::Cursor, fields, ::Type{TFriendElement_2})
fields[:fTreeName] = readtype(c.io, String)
fields[:fOwnFile] = readtype(c.io, Bool)
end

abstract type TAttBox2D <: ROOTStreamedObject end
struct TAttBox2D_0 <: TAttBox2D end
readfields!(c::Cursor, fields, ::Type{TAttBox2D_0}) = nothing

abstract type TText <: ROOTStreamedObject end
Base.@kwdef struct TText_3 <: TText
# TNamed
# TAttText
# TAttBox2D
fX::Float64
fY::Float64
end
function readfields!(c::Cursor, fields, ::Type{TText_3})
stream!(c, fields, TNamed)
stream!(c, fields, TAttText)
stream!(c, fields, TAttBox2D)
fields[:fX] = readtype(c.io, Float64)
fields[:fY] = readtype(c.io, Float64)
end

abstract type TLatex <: ROOTStreamedObject end
Base.@kwdef struct TLatex_2 <: TLatex
cursor::Cursor
# TText
fName::String
fTitle::String
fTextAngle::Float32
fTextSize::Float32
fTextAlign::Int16
fTextColor::Int16
fTextFont::Int16
fX::Float64
fY::Float64

# TAttLine
fLineColor::Int16
fLineStyle::Int16
fLineWidth::Int16

fLimitFactorSize::Int32
fOriginSize::Float64
end
function readfields!(c::Cursor, fields, ::Type{TLatex_2})
stream!(c, fields, TText)
stream!(c.io, fields, TAttLine) # TODO: define all these methods for the Cursor!
fields[:fLimitFactorSize] = readtype(c.io, Int32)
fields[:fOriginSize] = readtype(c.io, Float64)
end

abstract type TAttText <: ROOTStreamedObject end
Base.@kwdef struct TAttText_2 <: TAttText
cursor::Cursor
fTextAngle::Float32
fTextSize::Float32
fTextAlign::Int16
fTextColor::Int16
fTextFont::Int16
end
function readfields!(c::Cursor, fields, ::Type{TAttText_2})
fields[:fTextAngle] = readtype(c.io, Float32)
fields[:fTextSize] = readtype(c.io, Float32)
fields[:fTextAlign] = readtype(c.io, Int16)
fields[:fTextColor] = readtype(c.io, Int16)
fields[:fTextFont] = readtype(c.io, Int16)
end

# abstract type TBox <: ROOTStreamedObject end
# struct TBox_3 <: TBox
# fields::Dict{Symbol, Any}
# end
# function readfields!(io::IO, fields, ::Type{TBox_3})
# # skiptobj(io)
# # skip(io, 1)
# # stream!(io, fields, TAttLine)
# # stream!(io, fields, TAttFill)
# # skip(io, 6) # TAttBBox2D with cnt (u4) and vers (u2)
# skip(io, 140)
# end

abstract type TPaveText <: ROOTStreamedObject end
Base.@kwdef struct TPaveText_2 <: TPaveText
cursor::Cursor
# TBox
fX1NDC::Float64
fY1NDC::Float64
fX2NDC::Float64
fY2NDC::Float64
fBorderSize::Int32
fInit::Int32
fShadowColor::Int32
fCornerRadius::Float64
fOption::String
fName::String

fLabel::String
fLongest::Int32
fMargin::Float32
fLines
end
function readfields!(c::Cursor, fields, ::Type{TPaveText_2})
io = c.io
tkey = c.tkey
refs = c.refs

skip(io, 82) # TODO: what is this??

# TPave
# stream!(io, fields, TBox) # below is TBox readout since we could not
# figure out how to read it (see above)
fields[:fX1NDC] = readtype(io, Float64)
fields[:fY1NDC] = readtype(io, Float64)
fields[:fX2NDC] = readtype(io, Float64)
fields[:fY2NDC] = readtype(io, Float64)
fields[:fBorderSize] = readtype(io, Int32)
fields[:fInit] = readtype(io, Int32)
fields[:fShadowColor] = readtype(io, Int32)
fields[:fCornerRadius] = readtype(io, Float64)
fields[:fOption] = readtype(io, String)
fields[:fName] = readtype(io, String)

stream!(c, fields, TAttText)

fields[:fLabel] = readtype(io, String)
fields[:fLongest] = readtype(io, Int32)
fields[:fMargin] = readtype(io, Float32)
fields[:fLines] = readobjany!(io, tkey, refs)
end

abstract type TVirtualPaveStats <: ROOTStreamedObject end
struct TVirtualPaveStats_0 <: TVirtualPaveStats end
readfields!(c::Cursor, fields, ::Type{TVirtualPaveStats_0}) = nothing

abstract type TPaveStats <: ROOTStreamedObject end
Base.@kwdef struct TPaveStats_5 <: TPaveStats
cursor::Cursor
# TPaveText
fX1NDC::Float64
fY1NDC::Float64
fX2NDC::Float64
fY2NDC::Float64
fBorderSize::Int32
fInit::Int32
fShadowColor::Int32
fCornerRadius::Float64
fOption::String
fName::String
fLabel::String
fLongest::Int32
fMargin::Float32
fLines

fTextAngle::Float32
fTextSize::Float32
fTextAlign::Int16
fTextColor::Int16
fTextFont::Int16

fOptFit::Int32
fOptStat::Int32
fFitFormat::String
fStatFormat::String
fParent
end
function readfields!(c::Cursor, fields, ::Type{TPaveStats_5})
tkey = c.tkey
refs = c.refs

stream!(c, fields, TPaveText)
stream!(c, fields, TVirtualPaveStats)
fields[:fOptFit] = readtype(c.io, Int32)
fields[:fOptStat] = readtype(c.io, Int32)
fields[:fFitFormat] = readtype(c.io, String)
fields[:fStatFormat] = readtype(c.io, String)
fields[:fParent] = readobjany!(c.io, tkey, refs)
end
5 changes: 4 additions & 1 deletion src/streamers.jl
Original file line number Diff line number Diff line change
Expand Up @@ -219,7 +219,10 @@ function readobjany!(io, tkey::TKey, refs)
if tag == 0
return missing
elseif tag == 1
error("Returning parent is not implemented yet")
# TODO: Returning parent is not implemented yet but
# "missing" works fine so far ;) Create an issue if
# you have to deal with a "missing" from here.
return missing
elseif !haskey(refs, tag)
# skipping
seek(io, origin(tkey) + beg + bcnt + 4)
Expand Down
38 changes: 38 additions & 0 deletions test/histograms.jl
Original file line number Diff line number Diff line change
Expand Up @@ -78,4 +78,42 @@ using FHist
binlabels = ["Root", "Weight", "Preselection", "SelectGenPart", "GoodRunsList", "EventFilters", "SelectLeptons", "SelectJets", "Trigger", "ObjectsSelection", "SSPreselection", "NjetGeq4", "AK4CategTagHiggsJets", "AK4CategTagVBSJets", "AK4CategChannels", "AK4CategPresel"]
@test f["AK4CategPresel_cutflow"][:fXaxis_fModLabs].objects == binlabels
close(f)

f = UnROOT.samplefile("TH3F.root")
h = f["histogram"]
@test 45699.0 h[:fTsumw]
@test 45699.0 h[:fTsumw2]
@test 133.84083137224275 h[:fTsumwx]
@test 15255.146094082058 h[:fTsumwx2]
@test 489.6839597196363 h[:fTsumwy]
@test 60853.5600420268 h[:fTsumwy2]
@test 206.84178180322968 h[:fTsumwxy]
@test -424.5948793991111 h[:fTsumwz]
@test 137144.65414345643 h[:fTsumwz2]
@test -233.101696139249 h[:fTsumwxz]
@test -383.8438776730445 h[:fTsumwyz]
@test 2 == h[:fStatOverflows]
@test -1111 == h[:fMinimum]
@test -1111 == h[:fMaximum]
@test 99999.0 == h[:fEntries]
@test 0.0 == h[:fNormFactor]

@test 10 == h[:fXaxis_fNbins]
@test -1.0 == h[:fXaxis_fXmin]
@test 1.0 == h[:fXaxis_fXmax]
@test 0 == h[:fXaxis_fFirst]
@test 0 == h[:fXaxis_fLast]

@test 20 == h[:fYaxis_fNbins]
@test -2.0 == h[:fYaxis_fXmin]
@test 2.0 == h[:fYaxis_fXmax]
@test 0 == h[:fYaxis_fFirst]
@test 0 == h[:fYaxis_fLast]

@test 30 == h[:fZaxis_fNbins]
@test -3.0 == h[:fZaxis_fXmin]
@test 3.0 == h[:fZaxis_fXmax]
@test 0 == h[:fZaxis_fFirst]
@test 0 == h[:fZaxis_fLast]
close(f)
end
28 changes: 28 additions & 0 deletions test/samples/TH3F.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,28 @@
import random
import ROOT

# Create a ROOT file to store histograms
file = ROOT.TFile("TH3F.root", "RECREATE")

# Create a 3D histogram
histogram = ROOT.TH3F(
"histogram",
"Example 3D Histogram;X Axis;Y Axis;Z Axis",
10, -1, 1,
20, -2, 2,
30, -3, 3)

random.seed(42)

# Fill the histogram with some random data
for i in range(1, 100000):
x = random.random() * 3 - 1.5
y = random.random() * 5 - 2.5
z = random.random() * 7 - 3.5
histogram.Fill(x, y, z)

# Write the histogram to the ROOT file
histogram.Write()

# Close the ROOT file
file.Close()
Binary file added test/samples/TH3F.root
Binary file not shown.

0 comments on commit 3aad70a

Please sign in to comment.