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

Reading of c-arrays + more customstructs #167

Closed
wants to merge 5 commits into from
Closed
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions src/UnROOT.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,7 @@ module UnROOT
using LazyArrays
import Mmap: mmap
export ROOTFile, LazyBranch, LazyTree
#export FixSizeArray, FixSizeMatrix, FixSizeVector # would make it more readable in printed tables but not necessary

import Base: close, keys, get, getindex, getproperty, show, length, iterate, position, ntoh, reinterpret
ntoh(b::Bool) = b
Expand Down
38 changes: 24 additions & 14 deletions src/custom.jl
Original file line number Diff line number Diff line change
Expand Up @@ -31,23 +31,33 @@ end
# Custom struct interpretation
abstract type CustomROOTStruct end

struct FixLenVector{N, T} <: AbstractVector{T}
vec::SVector{N, T}
end
(::Type{FixLenVector{N, T}})() where {N, T} = FixLenVector(zero(SVector{N, T}))
Base.length(x::FixLenVector) = length(x.vec)
Base.length(::Type{FixLenVector{N, T}}) where {N, T} = N
Base.size(x::FixLenVector) = size(x.vec)
Base.eltype(x::FixLenVector) = eltype(x.vec)
Base.iterate(x::FixLenVector) = iterate(x.vec)
Base.iterate(x::FixLenVector, n) = iterate(x.vec, n)
Base.getindex(x::FixLenVector, n) = getindex(x.vec, n)
function Base.reinterpret(::Type{FixLenVector{N, T}}, v::AbstractVector{UInt8}) where {N, T}
# Fixed-size c-style array
struct FixSizeArray{S <: Tuple, T, N} <: AbstractArray{T, N}
arr::SArray{S, T, N}
end
const FixSizeVector{N, T} = FixSizeArray{Tuple{N}, T, 1}
const FixSizeMatrix{M, N, T} = FixSizeArray{Tuple{M, N}, T, 2}

(::Type{FixSizeArray{S, T, N}})() where {S, T, N} = FixSizeArray(zero(SArray{S, T, N}))
Base.length(x::FixSizeArray) = length(x.arr)
Base.length(t::Type{<:FixSizeArray}) = length(fieldtype(t,:arr))
Base.size(x::FixSizeArray) = size(x.arr)
Base.eltype(x::FixSizeArray) = eltype(x.arr)
Base.iterate(x::FixSizeArray) = iterate(x.arr)
Base.iterate(x::FixSizeArray, n) = iterate(x.arr, n)
Base.getindex(x::FixSizeArray, n...) = getindex(x.arr, n...)
# need this workaround since Base.sizeof is only correctly working for SArray
# with all 4 template arguments (note only 3 are used in FixSizeArray struct).
# Could get rid off this if it would be possible to defer the total size in our
# wrapper struct i.e. see https://github.com/JuliaLang/julia/issues/18466
Base.sizeof(::Type{<:FixSizeArray{S,T,N}}) where {S, T, N} = sizeof(T)*prod(fieldtypes(S))
function Base.reinterpret(::Type{<:FixSizeArray{S,T,N}}, v::AbstractVector{UInt8}) where {S, T, N}
vs = reinterpret(T, v)
@. vs = ntoh(vs)
FixLenVector(SVector{N, T}(vs))
FixSizeArray(SArray{S, T, N}(vs))
end
function interped_data(rawdata, rawoffsets, ::Type{T}, ::Type{Nojagg}) where {T <: FixLenVector}

function interped_data(rawdata, rawoffsets, ::Type{T}, ::Type{Nojagg}) where {T <: FixSizeArray}
n = sizeof(T)
[
reinterpret(T, x) for x in Base.Iterators.partition(rawdata, n)
Expand Down
8 changes: 4 additions & 4 deletions src/iteration.jl
Original file line number Diff line number Diff line change
Expand Up @@ -23,8 +23,8 @@ end

function array(f::ROOTFile, branch; raw=false)
ismissing(branch) && error("No branch found at $path")
(!raw && length(branch.fLeaves.elements) > 1) && error(
"Branches with multiple leaves are not supported yet. Try reading with `array(...; raw=true)`.",
!raw && (!haskey(f.customstructs, branch.fName) && length(branch.fLeaves.elements) > 1) && error(
"Branches with multiple leaves are not supported yet. Try reading with `array(...; raw=true)` or by including this branch into `customstructs`.",
)

rawdata, rawoffsets = readbranchraw(f, branch)
Expand All @@ -49,8 +49,8 @@ function basketarray(f::ROOTFile, path::AbstractString, ithbasket)
end
function basketarray(f::ROOTFile, branch, ithbasket)
ismissing(branch) && error("No branch found at $path")
length(branch.fLeaves.elements) > 1 && error(
"Branches with multiple leaves are not supported yet. Try reading with `array(...; raw=true)`.",
(!haskey(f.customstructs, branch.fName) && length(branch.fLeaves.elements) > 1) && error(
"Branches with multiple leaves are not supported yet. Try reading with `array(...; raw=true)` or by including this branch into `customstructs`.",
)

rawdata, rawoffsets = readbasket(f, branch, ithbasket)
Expand Down
9 changes: 5 additions & 4 deletions src/root.jl
Original file line number Diff line number Diff line change
Expand Up @@ -318,9 +318,9 @@ function auto_T_JaggT(f::ROOTFile, branch; customstructs::Dict{String, Type})
leaf = first(branch.fLeaves.elements)
_type = Nothing
_jaggtype = JaggType(f, branch, leaf)
if hasproperty(branch, :fClassName)
classname = branch.fClassName # the C++ class name, such as "vector<int>"
parentname = branch.fParentName # assuming it has a parent ;)
if hasproperty(branch, :fClassName) || haskey(customstructs, branch.fName)
classname = hasproperty(branch, :fClassName) ? branch.fClassName : branch.fName # the C++ class name, such as "vector<int>" if C++-class
parentname = hasproperty(branch, :fParentName) ? branch.fParentName : nothing # assuming it has a parent ;)
try
# this will call a customize routine if defined by user
# see custom.jl
Expand Down Expand Up @@ -394,7 +394,8 @@ function auto_T_JaggT(f::ROOTFile, branch; customstructs::Dict{String, Type})
else
_type = primitivetype(leaf)
if leaf.fLen > 1 # treat NTuple as Nojagg since size is static
_type = FixLenVector{Int(leaf.fLen), _type}
dims = LeafDims(leaf)
_type = FixSizeArray{Tuple{Int.(dims)...}, _type, length(dims)}
return _type, Nojagg
end
_type = _jaggtype === Nojagg ? _type : Vector{_type}
Expand Down
8 changes: 8 additions & 0 deletions src/utils.jl
Original file line number Diff line number Diff line change
Expand Up @@ -104,3 +104,11 @@ end
function samplefile(filename::AbstractString)
return ROOTFile(normpath(joinpath(@__DIR__, "../test/samples", filename)))
end

function LeafDims(leaf)
try
parse.(Int,split(strip(match(r"\[.*\]",leaf.fTitle).match,['[',']']),"]["))
catch
leaf.fLen
end
end