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

TATTMarker_3 not defined in UnROOT #375

Closed
pranphy opened this issue Jan 16, 2025 · 5 comments · Fixed by #376
Closed

TATTMarker_3 not defined in UnROOT #375

pranphy opened this issue Jan 16, 2025 · 5 comments · Fixed by #376

Comments

@pranphy
Copy link

pranphy commented Jan 16, 2025

I have a root file with multiple root histograms. I am trying to read it with unroot but its failing with error TAttMarker_3 not defined in UnROOT.

rootfile = "sync-edge-check6.root";
f = ROOTFile(rootfile)
sync-edge-check6.root
└─ synch_check (TDirectory)
   └─ step_field (TDirectory)
      ├─ syncyed (TH1D)
      ├─ syncxyy (TH2D)
      ├─ syncxyd1 (TH2D)
      ├─ syncvz (TH2D)
      ├─ syncvzd1 (TH2D)
      └─ synce (TH1D)

When I try to access the histogram like this

f["synch_check/step_field/syncyed"]
ERROR: UndefVarError: `TAttMarker_3` not defined in `UnROOT`
Stacktrace:
  [1] stream!(io::IOBuffer, fields::Dict{Symbol, Any}, ::Type{UnROOT.TAttMarker}; check::Bool)
    @ UnROOT ~/.julia/packages/UnROOT/ABV18/src/streamers.jl:590
  [2] stream!
    @ ~/.julia/packages/UnROOT/ABV18/src/streamers.jl:586 [inlined]
  [3] TH(io::UnROOT.MmapStream, tkey::UnROOT.TKey32, refs::Dict{Int32, Any})
    @ UnROOT ~/.julia/packages/UnROOT/ABV18/src/bootstrap.jl:984
  [4] TH1D(io::UnROOT.MmapStream, tkey::UnROOT.TKey32, refs::Dict{Int32, Any})
    @ UnROOT ~/.julia/packages/UnROOT/ABV18/src/bootstrap.jl:954
  [5] getindex(d::UnROOT.ROOTDirectory, s::String)
    @ UnROOT ~/.julia/packages/UnROOT/ABV18/src/root.jl:195
  [6] getindex(d::UnROOT.ROOTDirectory, s::String)
    @ UnROOT ~/.julia/packages/UnROOT/ABV18/src/root.jl:191
  [7] _getindex(f::ROOTFile, s::String)
    @ UnROOT ~/.julia/packages/UnROOT/ABV18/src/root.jl:171
  [8] #150
    @ ~/.julia/packages/UnROOT/ABV18/src/root.jl:163 [inlined]
  [9] get!(default::UnROOT.var"#150#151"{ROOTFile, String}, h::Dict{Any, Any}, key::String)
    @ Base ./dict.jl:458
 [10] getindex(f::ROOTFile, s::String)
    @ UnROOT ~/.julia/packages/UnROOT/ABV18/src/root.jl:162
 [11] top-level scope
    @ none:1

The file: sync-edge-check6.root

Edit: After some investigation I found that it happened only to my files that I created recently with root. Turns this issue only happens for the root files that I created with my latest upgraded root "6.34" I downgraded root back to 6.32.04 and recreated this rootfile and UnROOT correctly parsed the histogram.

@Moelf
Copy link
Member

Moelf commented Jan 16, 2025

latest upgraded root "6.34" I downgraded root back to 6.32.04 and recreated this rootfile and UnROOT correctly parsed the histogram.

when this happens it means ROOT has internally changed / updated the C++ side definition of TH1D

sigh, I thought they would leave such an old object alone but we should be able to fix this. cc. @tamasgal

@tamasgal
Copy link
Member

Ah OK, so this should be it I think: root-project/root@ca1f047

I'll have a look.

@tamasgal
Copy link
Member

@pranphy sorry for the long delay, version 0.10.34 is on the way with the fix!

julia> using UnROOT

julia> f = ROOTFile("/Users/tamasgal/Downloads/sync-edge-check6.root")
ROOTFile with 1 entry and 14 streamers.
/Users/tamasgal/Downloads/sync-edge-check6.root
└─ synch_check (TDirectory)
   └─ step_field (TDirectory)
      ├─ syncyed (TH1D)
      ├─ syncxyy (TH2D)
      ├─ syncxyd1 (TH2D)
      ├─ syncvz (TH2D)
      ├─ syncvzd1 (TH2D)
      └─ synce (TH1D)


julia> f["synch_check/step_field/syncyed"]
Dict{Symbol, Any} with 102 entries:
  :fZaxis_fTitleSize   => 0.035
  :fYaxis_fXmin        => 0.0
  :fYaxis_fNbins       => 1
  :fZaxis_fLast        => 0
  :fYaxis_fName        => "yaxis"
  :fMaximum            => -1111.0
  :fZaxis_fLabelFont   => 42
  :fZaxis_fLabelSize   => 0.035
  :fZaxis_fXmin        => 0.0
  :fZaxis_fModLabs     => missing
  :fYaxis_fLabelColor  => 1
  :fZaxis_fLabels      => missing
  :fYaxis_fFirst       => 0
  :fYaxis_fTitleSize   => 0.035
  :fLineColor          => 602
  :fFillColor          => 0
  :fYaxis_fTimeFormat  => ""
  :fNcells             => 102
  :fYaxis_fXbins       => Float64[]
  :fXaxis_fLabelOffset => 0.005
  :fZaxis_fXbins       => Float64[]
  :fZaxis_fName        => "zaxis"
  :fTitle              => "Synchrotron Photon"
  :fTsumw2             => 15545.0
  :fXaxis_fTickLength  => 0.03
  :fXaxis_fTimeFormat  => ""
  :fZaxis_fLabelOffset => 0.005
  :fMarkerStyle        => 1
  :fContour            => Float64[]
  :fZaxis_fXmax        => 1.0
  :fZaxis_fTitleColor  => 1
  :fBarOffset          => 0
  :fXaxis_fTitleOffset => 1.0
  :fXaxis_fModLabs     => missing
  :fYaxis_fLabels      => missing
  :fXaxis_fTitleSize   => 0.035
  :fYaxis_fNdivisions  => 510
  :fXaxis_fTitleColor  => 1
  :fZaxis_fAxisColor   => 1
  :fBarWidth           => 1000
  :fYaxis_fTimeDisplay => false
  :fLineWidth          => 1
  :fMarkerColor        => 1
  :fYaxis_fTitle       => " Count"
  :fXaxis_fNbins       => 100
  :fTsumwx2            => 5.96962e8
  :fXaxis_fNdivisions  => 510
  :fYaxis_fLabelOffset => 0.005
  :fZaxis_fTimeFormat  => ""
                      => 

julia>

@tamasgal
Copy link
Member

For the sake of completeness, here is how to get the histogram data ;)

julia> h = UnROOT.parseTH(f["synch_check/step_field/syncyed"]; raw=false)
edges: [178.0, 178.42, 178.84, 179.26, 179.68, 180.1, 180.52, 180.94, 181.36, 181.78    216.22, 216.64, 217.06, 217.48, 217.9, 218.32, 218.74, 219.16, 219.58, 220.0]
bin counts: [262.0, 230.0, 223.0, 180.0, 209.0, 197.0, 190.0, 178.0, 197.0, 208.0    33.0, 27.0, 14.0, 6.0, 9.0, 6.0, 2.0, 4.0, 3.0, 0.0]
total count: 15549.0

julia> typeof(h)
FHist.Hist1D{Float64}

julia> h = UnROOT.parseTH(f["synch_check/step_field/syncvzd1"]; raw=false)
edges: ([-8000.0, -7840.0, -7680.0, -7520.0, -7360.0, -7200.0, -7040.0, -6880.0, -6720.0, -6560.0    6560.0, 6720.0, 6880.0, 7040.0, 7200.0, 7360.0, 7520.0, 7680.0, 7840.0, 8000.0], [0.0, 3.0, 6.0, 9.0, 12.0, 15.0, 18.0, 21.0, 24.0, 27.0    273.0, 276.0, 279.0, 282.0, 285.0, 288.0, 291.0, 294.0, 297.0, 300.0])
bin counts: [0.0 0.0  0.0 0.0; 0.0 0.0  0.0 0.0;  ; 0.0 0.0  0.0 0.0; 0.0 0.0  0.0 0.0]
total count: 16430.0

julia> typeof(h)
FHist.Hist2D{Float64}

@pranphy
Copy link
Author

pranphy commented Jan 23, 2025

@tamasgal Thank you very much. Really appreciate this.

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

Successfully merging a pull request may close this issue.

3 participants