Skip to content

Commit 8b16c8f

Browse files
authored
handle [var][fix] style doubly jagged branch (#188)
* handle [var][fix] style doubly jagged branch * fix tests for 1.6
1 parent 232e85b commit 8b16c8f

File tree

6 files changed

+99
-12
lines changed

6 files changed

+99
-12
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,7 +1,7 @@
11
name = "UnROOT"
22
uuid = "3cd96dde-e98d-4713-81e9-a4a1b0235ce9"
33
authors = ["Tamas Gal", "Jerry Ling", "Johannes Schumann", "Nick Amin"]
4-
version = "0.8.15"
4+
version = "0.8.16"
55

66
[deps]
77
AbstractTrees = "1520ce14-60c1-5f80-bbc7-55ef81b5835c"

src/custom.jl

Lines changed: 8 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -43,8 +43,16 @@ Base.iterate(x::FixLenVector) = iterate(x.vec)
4343
Base.iterate(x::FixLenVector, n) = iterate(x.vec, n)
4444
Base.getindex(x::FixLenVector, n) = getindex(x.vec, n)
4545

46+
#N.B this is a hack, we deal with ntoh in reinterpret step
47+
Base.ntoh(x::FixLenVector) = x
4648
const VorView = Union{Vector{UInt8}, SubArray{UInt8, 1, Vector{UInt8}, Tuple{UnitRange{Int64}}, true}}
4749

50+
function Base.reinterpret(::Type{Vector{FixLenVector{N, T}}}, data::Vector{UInt8}) where {N,T}
51+
vs = reinterpret(T, data)
52+
@. vs = ntoh(vs)
53+
return FixLenVector.(reinterpret(SVector{N, T}, vs))
54+
end
55+
4856
function Base.reinterpret(::Type{FixLenVector{N, T}}, v::VorView) where {N, T}
4957
vs = reinterpret(T, v)
5058
@. vs = ntoh(vs)

src/root.jl

Lines changed: 22 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -407,14 +407,31 @@ function auto_T_JaggT(f::ROOTFile, branch; customstructs::Dict{String, Type})
407407
end
408408
end
409409
else
410+
# since no classname were found, we now try to determine
411+
# type based on leaf information
412+
_type, _jaggtype = leaf_jaggtype(leaf, _jaggtype)
413+
end
414+
return _type, _jaggtype
415+
end
416+
417+
function leaf_jaggtype(leaf, _jaggtype)
410418
_type = primitivetype(leaf)
411-
if leaf.fLen > 1 # treat NTuple as Nojagg since size is static
412-
_type = FixLenVector{Int(leaf.fLen), _type}
413-
return _type, Nojagg
419+
leafLen = leaf.fLen
420+
if leafLen > 1 # treat NTuple as Nojagg since size is static
421+
_fTitle = replace(leaf.fTitle, "[$(leafLen)]" => "")
422+
# looking for more `[var]`
423+
m = match(r"\[\D+\]", _fTitle)
424+
_type = FixLenVector{Int(leafLen), _type}
425+
if isnothing(m)
426+
return _type, Nojagg
427+
else
428+
#FIXME this only handles [var][fix] case
429+
return Vector{_type}, Nooffsetjagg
430+
end
414431
end
415432
_type = _jaggtype === Nojagg ? _type : Vector{_type}
416-
end
417-
return _type, _jaggtype
433+
434+
return _type, _jaggtype
418435
end
419436

420437

test/runtests.jl

Lines changed: 41 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -318,8 +318,9 @@ end
318318
@test v2[2].int32_array == data[3].int32_array
319319
end
320320

321-
@testset "Doubly jagged branches" begin
322-
rootfile = ROOTFile(joinpath(SAMPLES_DIR, "tree_with_doubly_jagged.root"))
321+
@testset "Doubly jagged [var][var] branches" begin
322+
# this is vector<vector<blah>>
323+
rootfile = UnROOT.samplefile("tree_with_doubly_jagged.root")
323324
vvi = [[[2], [3, 5]], [[7, 9, 11], [13]], [[17], [19], []], [], [[]]]
324325
vvf = [[[2.5], [3.5, 5.5]], [[7.5, 9.5, 11.5], [13.5]], [[17.5], [19.5], []], [], [[]]]
325326
@test UnROOT.array(rootfile, "t1/bi") == vvi
@@ -331,6 +332,33 @@ end
331332
close(rootfile)
332333
end
333334

335+
@testset "Doubly jagged [var][fix] branches" begin
336+
# issue #187
337+
# this is vector<Int[N]>
338+
f = UnROOT.samplefile("tree_with_varfix_doubly_jagged.root")
339+
tree = LazyTree(f, "outtree")
340+
@test tree.nparticles == [4,3,2]
341+
@test length.(tree.P) == [4,3,2]
342+
@test eltype(tree.P[1]) <: AbstractVector
343+
# also compared to uproot
344+
@test tree[1].P == [
345+
[0.9411764705882353, 0.8888888888888888, 0.8421052631578947, 0.8],
346+
[1.0, 0.9285714285714286, 0.8666666666666667, 0.8125],
347+
[1.1111111111111112, 1.0, 0.9090909090909091, 0.8333333333333334],
348+
[1.4, 1.1666666666666667, 1.0, 0.875]
349+
]
350+
@test tree[3].P == [
351+
[0.8222222222222222,
352+
0.8043478260869565,
353+
0.7872340425531915,
354+
0.7708333333333334],
355+
[0.8292682926829268,
356+
0.8095238095238095,
357+
0.7906976744186046,
358+
0.7727272727272727]
359+
]
360+
end
361+
334362
@testset "NanoAOD" begin
335363
rootfile = ROOTFile(joinpath(SAMPLES_DIR, "NanoAODv5_sample.root"))
336364
event = UnROOT.array(rootfile, "Events/event")
@@ -641,8 +669,11 @@ nthreads == 1 && @warn "Running on a single thread. Please re-run the test suite
641669

642670

643671
if get(ENV, "CI", "false") == "true"
644-
# Make sure CI runs with more than 1 thread
645-
@test Threads.nthreads()>1 skip=(nthreads == 1)
672+
if nthreads >= 1
673+
@test Threads.nthreads()>1
674+
else
675+
@warn "CI wasn't run with multi thread"
676+
end
646677
end
647678
nmus = zeros(Int, Threads.nthreads())
648679
Threads.@threads for i in 1:length(t)
@@ -671,15 +702,19 @@ t = LazyTree(ROOTFile(joinpath(SAMPLES_DIR, "NanoAODv5_sample.root")), "Events",
671702
Threads.@threads for evt in t
672703
nmus[Threads.threadid()] += length(evt.Muon_pt)
673704
end
674-
@test count(>(0), nmus) > 1 skip = (nthreads==1)# test @threads is actually threading
705+
if nthreads > 1
706+
@test count(>(0), nmus) > 1# test @threads is actually threading
707+
end
675708
@test sum(nmus) == 878
676709

677710

678711
nmus .= 0
679712
Threads.@threads for evt in t
680713
nmus[Threads.threadid()] += length(evt.Muon_pt)
681714
end
682-
@test count(>(0), nmus) > 1 skip = (nthreads==1)# test @threads is actually threading
715+
if nthreads > 1
716+
@test count(>(0), nmus) > 1
717+
end
683718
@test sum(nmus) == 878
684719

685720
nmus .= 0
Lines changed: 27 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,27 @@
1+
#include "TFile.h"
2+
#include "TTree.h"
3+
4+
int maketree(){
5+
TFile f("tree_with_varfix_doubly_jagged.root", "RECREATE", "");
6+
TTree tree = TTree("outtree", "outtree");
7+
int nparticles{};
8+
double P[100][4];
9+
tree.Branch("nparticles", &nparticles, "nparticles/I");
10+
tree.Branch("P", P, "P[nparticles][4]/D");
11+
double counter1 = 1;
12+
double counter2 = 1;
13+
for (auto ev = 0; ev<3; ++ev){
14+
nparticles = 4-ev;
15+
for (auto i = nparticles; i>=0; --i){
16+
counter1 += 3;
17+
for (auto j = 0; j<=3; ++j){
18+
P[i][j] = counter1 / (counter2);
19+
counter2++;
20+
}
21+
}
22+
tree.Fill();
23+
}
24+
f.Write();
25+
f.Close();
26+
return 0;
27+
}
5.92 KB
Binary file not shown.

0 commit comments

Comments
 (0)