EDM4hep.AnalysisModule

Analysis module for EDM4hep.jl

A number of data structures and functions to support analysis of EDM4hep data in multithreaded mode.

Base.append!Method
Base.append!(d1::AbstractAnalysisData, d2::AbstractAnalysisData)

Default function to reset the user analysis data structure in case is not provided explicitely.

Base.empty!Method
Base.empty!(data::AbstractAnalysisData)

Default function to reset the user analysis data structure in case is not provided explicitely.

EDM4hep.Analysis.do_analysis!Method

doanalysis!(data::AbstractAnalysisData, analysis, reader, events; mt::Bool=false, tasksper_thread::Int=4)

Perform an analysis on all events by executing the analysis function. The iteration will be chunked and distributed to diffrent tasks running on different threads if the option argument mt is set to true. The results in data for all the cbhunks will be merged at the end of the analysis.

EDM4hep.EDM4hepModule

Main module for EDM4hep.jl – Key4hep Event Data Model for Julia.

All data model types are exported from this module for public use

Exports

EDM4hep.CaloHitContributionType

Monte Carlo contribution to SimCalorimeterHit

  • Author: F.Gaede, DESY

Fields

  • PDG::Int32: PDG code of the shower particle that caused this contribution.
  • energy::Float32: energy in [GeV] of the this contribution
  • time::Float32: time in [ns] of this contribution
  • stepPosition::Vector3f: position of this energy deposition (step) [mm]

Relations

  • particle::MCParticle: primary MCParticle that caused the shower responsible for this contribution to the hit.
EDM4hep.CalorimeterHitType

Calorimeter hit

  • Author: F.Gaede, DESY

Fields

  • cellID::UInt64: detector specific (geometrical) cell id.
  • energy::Float32: energy of the hit in [GeV].
  • energyError::Float32: error of the hit energy in [GeV].
  • time::Float32: time of the hit in [ns].
  • position::Vector3f: position of the hit in world coordinates in [mm].
  • type::Int32: type of hit. Mapping of integer types to names via collection parameters "CalorimeterHitTypeNames" and "CalorimeterHitTypeValues".
EDM4hep.ClusterType

Calorimeter Hit Cluster

  • Author: F.Gaede, DESY

Fields

  • type::Int32: flagword that defines the type of cluster. Bits 16-31 are used internally.
  • energy::Float32: energy of the cluster [GeV]
  • energyError::Float32: error on the energy
  • position::Vector3f: position of the cluster [mm]
  • positionError::SVector{6,Float32}: covariance matrix of the position (6 Parameters)
  • iTheta::Float32: intrinsic direction of cluster at position Theta. Not to be confused with direction cluster is seen from IP.
  • phi::Float32: intrinsic direction of cluster at position - Phi. Not to be confused with direction cluster is seen from IP.
  • directionError::Vector3f: covariance matrix of the direction (3 Parameters) [mm^2]
  • shapeParameters::PVector{Float32}: shape parameters - check/set collection parameter ClusterShapeParameters for size and names of parameters.
  • subdetectorEnergies::PVector{Float32}: energy observed in a particular subdetector. Check/set collection parameter ClusterSubdetectorNames for decoding the indices of the array.

Relations

  • clusters::Cluster: clusters that have been combined to this cluster.
  • hits::CalorimeterHit: hits that have been combined to this cluster.
  • particleIDs::ParticleID: particle IDs (sorted by their likelihood)

Methods

  • setShapeParameters(object::Cluster, v::AbstractVector{Float32}): assign a set of values to the shapeParameters vector member
  • setSubdetectorEnergies(object::Cluster, v::AbstractVector{Float32}): assign a set of values to the subdetectorEnergies vector member
  • pushToClusters(obj::Cluster, robj::Cluster): push related object to the clusters relation
  • popFromClusters(obj::Cluster): pop last related object from clusters relation
  • pushToHits(obj::Cluster, robj::CalorimeterHit): push related object to the hits relation
  • popFromHits(obj::Cluster): pop last related object from hits relation
  • pushToParticleIDs(obj::Cluster, robj::ParticleID): push related object to the particleIDs relation
  • popFromParticleIDs(obj::Cluster): pop last related object from particleIDs relation
EDM4hep.EventHeaderType

Event Header. Additional parameters are assumed to go into the metadata tree.

  • Author: F.Gaede

Fields

  • eventNumber::Int32: event number
  • runNumber::Int32: run number
  • timeStamp::UInt64: time stamp
  • weight::Float32: event weight
EDM4hep.HitLevelDataType

HitLevelData

Fields

  • cellID::UInt64: cell id
  • N::UInt32: number of reconstructed ionization cluster.
  • eDep::Float32: reconstructed energy deposit [GeV].
  • pathLength::Float32: track path length [mm].
EDM4hep.HypothesisType

Hypothesis

Fields

  • chi2::Float32: chi2
  • expected::Float32: expected value
  • sigma::Float32: sigma value
EDM4hep.MCParticleType

The Monte Carlo particle - based on the lcio::MCParticle.

  • Author: F.Gaede, DESY

Fields

  • PDG::Int32: PDG code of the particle
  • generatorStatus::Int32: status of the particle as defined by the generator
  • simulatorStatus::Int32: status of the particle from the simulation program - use BIT constants below
  • charge::Float32: particle charge
  • time::Float32: creation time of the particle in [ns] wrt. the event, e.g. for preassigned decays or decays in flight from the simulator.
  • mass::Float64: mass of the particle in [GeV]
  • vertex::Vector3d: production vertex of the particle in [mm].
  • endpoint::Vector3d: endpoint of the particle in [mm]
  • momentum::Vector3d: particle 3-momentum at the production vertex in [GeV]
  • momentumAtEndpoint::Vector3d: particle 3-momentum at the endpoint in [GeV]
  • spin::Vector3f: spin (helicity) vector of the particle.
  • colorFlow::Vector2i: color flow as defined by the generator

Relations

  • parents::MCParticle: The parents of this particle.
  • daughters::MCParticle: The daughters this particle.

Methods

  • pushToParents(obj::MCParticle, robj::MCParticle): push related object to the parents relation
  • popFromParents(obj::MCParticle): pop last related object from parents relation
  • pushToDaughters(obj::MCParticle, robj::MCParticle): push related object to the daughters relation
  • popFromDaughters(obj::MCParticle): pop last related object from daughters relation
EDM4hep.MCRecoCaloAssociationType

Association between a CaloHit and the corresponding simulated CaloHit

  • Author: C. Bernet, B. Hegner

Fields

  • weight::Float32: weight of this association

Relations

  • rec::CalorimeterHit: reference to the reconstructed hit
  • sim::SimCalorimeterHit: reference to the simulated hit
EDM4hep.MCRecoCaloParticleAssociationType

Association between a CalorimeterHit and a MCParticle

  • Author: Placido Fernandez Declara

Fields

  • weight::Float32: weight of this association

Relations

  • rec::CalorimeterHit: reference to the reconstructed hit
  • sim::MCParticle: reference to the Monte-Carlo particle
EDM4hep.MCRecoClusterParticleAssociationType

Association between a Cluster and a MCParticle

  • Author: Placido Fernandez Declara

Fields

  • weight::Float32: weight of this association

Relations

  • rec::Cluster: reference to the cluster
  • sim::MCParticle: reference to the Monte-Carlo particle
EDM4hep.MCRecoParticleAssociationType

Used to keep track of the correspondence between MC and reconstructed particles

  • Author: C. Bernet, B. Hegner

Fields

  • weight::Float32: weight of this association

Relations

  • rec::ReconstructedParticle: reference to the reconstructed particle
  • sim::MCParticle: reference to the Monte-Carlo particle
EDM4hep.MCRecoTrackParticleAssociationType

Association between a Track and a MCParticle

  • Author: Placido Fernandez Declara

Fields

  • weight::Float32: weight of this association

Relations

  • rec::Track: reference to the track
  • sim::MCParticle: reference to the Monte-Carlo particle
EDM4hep.MCRecoTrackerAssociationType

Association between a TrackerHit and the corresponding simulated TrackerHit

  • Author: C. Bernet, B. Hegner

Fields

  • weight::Float32: weight of this association

Relations

  • rec::TrackerHit: reference to the reconstructed hit
  • sim::SimTrackerHit: reference to the simulated hit
EDM4hep.MCRecoTrackerHitPlaneAssociationType

Association between a TrackerHitPlane and the corresponding simulated TrackerHit

  • Author: Placido Fernandez Declara

Fields

  • weight::Float32: weight of this association

Relations

  • rec::TrackerHitPlane: reference to the reconstructed hit
  • sim::SimTrackerHit: reference to the simulated hit
EDM4hep.ParticleIDType

ParticleID

  • Author: F.Gaede, DESY

Fields

  • type::Int32: userdefined type
  • PDG::Int32: PDG code of this id - ( 999999 ) if unknown.
  • algorithmType::Int32: type of the algorithm/module that created this hypothesis
  • likelihood::Float32: likelihood of this hypothesis - in a user defined normalization.
  • parameters::PVector{Float32}: parameters associated with this hypothesis. Check/set collection parameters ParameterNames_PIDAlgorithmTypeName for decoding the indices.

Methods

  • setParameters(object::ParticleID, v::AbstractVector{Float32}): assign a set of values to the parameters vector member
EDM4hep.QuantityType

Quantity

Fields

  • type::Int16: flag identifying how to interpret the quantity
  • value::Float32: value of the quantity
  • error::Float32: error on the value of the quantity
EDM4hep.RawCalorimeterHitType

Raw calorimeter hit

  • Author: F.Gaede, DESY

Fields

  • cellID::UInt64: detector specific (geometrical) cell id.
  • amplitude::Int32: amplitude of the hit in ADC counts.
  • timeStamp::Int32: time stamp for the hit.
EDM4hep.RawTimeSeriesType

Raw data of a detector readout

  • Author: F.Gaede, DESY

Fields

  • cellID::UInt64: detector specific cell id.
  • quality::Int32: quality flag for the hit.
  • time::Float32: time of the hit [ns].
  • charge::Float32: integrated charge of the hit [fC].
  • interval::Float32: interval of each sampling [ns].
  • adcCounts::PVector{Int32}: raw data (32-bit) word at i.

Methods

  • setAdcCounts(object::RawTimeSeries, v::AbstractVector{Int32}): assign a set of values to the adcCounts vector member
EDM4hep.RecDqdxType

dN/dx or dE/dx info of Track.

  • Author: Wenxing Fang, IHEP

Fields

  • dQdx::Quantity: the reconstructed dEdx or dNdx and its error
  • particleType::Int16: particle type, e(0),mu(1),pi(2),K(3),p(4).
  • type::Int16: type.
  • hypotheses::SVector{5,Hypothesis}: 5 particle hypothesis
  • hitData::PVector{HitLevelData}: hit level data

Relations

  • track::Track: the corresponding track.

Methods

  • setHitData(object::RecDqdx, v::AbstractVector{HitLevelData}): assign a set of values to the hitData vector member
EDM4hep.RecIonizationClusterType

Reconstructed Ionization Cluster

  • Author: Wenxing Fang, IHEP

Fields

  • cellID::UInt64: cell id.
  • significance::Float32: significance.
  • type::Int16: type.

Relations

  • trackerPulse::TrackerPulse: the TrackerPulse used to create the ionization cluster.

Methods

  • pushToTrackerPulse(obj::RecIonizationCluster, robj::TrackerPulse): push related object to the trackerPulse relation
  • popFromTrackerPulse(obj::RecIonizationCluster): pop last related object from trackerPulse relation
EDM4hep.RecoParticleVertexAssociationType

Association between a Reconstructed Particle and a Vertex

  • Author: Placido Fernandez Declara

Fields

  • weight::Float32: weight of this association

Relations

  • rec::ReconstructedParticle: reference to the reconstructed particle
  • vertex::Vertex: reference to the vertex
EDM4hep.ReconstructedParticleType

Reconstructed Particle

  • Author: F.Gaede, DESY

Fields

  • type::Int32: type of reconstructed particle. Check/set collection parameters ReconstructedParticleTypeNames and ReconstructedParticleTypeValues.
  • energy::Float32: [GeV] energy of the reconstructed particle. Four momentum state is not kept consistent internally.
  • momentum::Vector3f: [GeV] particle momentum. Four momentum state is not kept consistent internally.
  • referencePoint::Vector3f: [mm] reference, i.e. where the particle has been measured
  • charge::Float32: charge of the reconstructed particle.
  • mass::Float32: [GeV] mass of the reconstructed particle, set independently from four vector. Four momentum state is not kept consistent internally.
  • goodnessOfPID::Float32: overall goodness of the PID on a scale of [0;1]
  • covMatrix::SVector{10,Float32}: cvariance matrix of the reconstructed particle 4vector (10 parameters). Stored as lower triangle matrix of the four momentum (px,py,pz,E), i.e. cov(px,px), cov(py,##

Relations

  • startVertex::Vertex: start vertex associated to this particle
  • particleIDUsed::ParticleID: particle Id used for the kinematics of this particle
  • clusters::Cluster: clusters that have been used for this particle.
  • tracks::Track: tracks that have been used for this particle.
  • particles::ReconstructedParticle: reconstructed particles that have been combined to this particle.
  • particleIDs::ParticleID: particle Ids (not sorted by their likelihood)

Methods

  • pushToClusters(obj::ReconstructedParticle, robj::Cluster): push related object to the clusters relation
  • popFromClusters(obj::ReconstructedParticle): pop last related object from clusters relation
  • pushToTracks(obj::ReconstructedParticle, robj::Track): push related object to the tracks relation
  • popFromTracks(obj::ReconstructedParticle): pop last related object from tracks relation
  • pushToParticles(obj::ReconstructedParticle, robj::ReconstructedParticle): push related object to the particles relation
  • popFromParticles(obj::ReconstructedParticle): pop last related object from particles relation
  • pushToParticleIDs(obj::ReconstructedParticle, robj::ParticleID): push related object to the particleIDs relation
  • popFromParticleIDs(obj::ReconstructedParticle): pop last related object from particleIDs relation
EDM4hep.SimCalorimeterHitType

Simulated calorimeter hit

  • Author: F.Gaede, DESY

Fields

  • cellID::UInt64: ID of the sensor that created this hit
  • energy::Float32: energy of the hit in [GeV].
  • position::Vector3f: position of the hit in world coordinates in [mm].

Relations

  • contributions::CaloHitContribution: Monte Carlo step contribution - parallel to particle

Methods

  • pushToContributions(obj::SimCalorimeterHit, robj::CaloHitContribution): push related object to the contributions relation
  • popFromContributions(obj::SimCalorimeterHit): pop last related object from contributions relation
EDM4hep.SimPrimaryIonizationClusterType

Simulated Primary Ionization

  • Author: Wenxing Fang, IHEP

Fields

  • cellID::UInt64: cell id.
  • time::Float32: the primary ionization's time in the lab frame [ns].
  • position::Vector3d: the primary ionization's position [mm].
  • type::Int16: type.
  • electronCellID::PVector{UInt64}: cell id.
  • electronTime::PVector{Float32}: the time in the lab frame [ns].
  • electronPosition::PVector{Vector3d}: the position in the lab frame [mm].
  • pulseTime::PVector{Float32}: the pulse's time in the lab frame [ns].
  • pulseAmplitude::PVector{Float32}: the pulse's amplitude [fC].

Relations

  • mcparticle::MCParticle: the particle that caused the ionizing collisions.

Methods

  • setElectronCellID(object::SimPrimaryIonizationCluster, v::AbstractVector{UInt64}): assign a set of values to the electronCellID vector member
  • setElectronTime(object::SimPrimaryIonizationCluster, v::AbstractVector{Float32}): assign a set of values to the electronTime vector member
  • setElectronPosition(object::SimPrimaryIonizationCluster, v::AbstractVector{Vector3d}): assign a set of values to the electronPosition vector member
  • setPulseTime(object::SimPrimaryIonizationCluster, v::AbstractVector{Float32}): assign a set of values to the pulseTime vector member
  • setPulseAmplitude(object::SimPrimaryIonizationCluster, v::AbstractVector{Float32}): assign a set of values to the pulseAmplitude vector member
EDM4hep.SimTrackerHitType

Simulated tracker hit

  • Author: F.Gaede, DESY

Fields

  • cellID::UInt64: ID of the sensor that created this hit
  • EDep::Float32: energy deposited in the hit [GeV].
  • time::Float32: proper time of the hit in the lab frame in [ns].
  • pathLength::Float32: path length of the particle in the sensitive material that resulted in this hit.
  • quality::Int32: quality bit flag.
  • position::Vector3d: the hit position in [mm].
  • momentum::Vector3f: the 3-momentum of the particle at the hits position in [GeV]

Relations

  • mcparticle::MCParticle: MCParticle that caused the hit.
EDM4hep.TimeSeriesType

Calibrated Detector Data

  • Author: Wenxing Fang, IHEP

Fields

  • cellID::UInt64: cell id.
  • time::Float32: begin time [ns].
  • interval::Float32: interval of each sampling [ns].
  • amplitude::PVector{Float32}: calibrated detector data.

Methods

  • setAmplitude(object::TimeSeries, v::AbstractVector{Float32}): assign a set of values to the amplitude vector member
EDM4hep.TrackType

Reconstructed track

  • Author: F.Gaede, DESY

Fields

  • type::Int32: flagword that defines the type of track.Bits 16-31 are used internally
  • chi2::Float32: Chi^2 of the track fit
  • ndf::Int32: number of degrees of freedom of the track fit
  • dEdx::Float32: dEdx of the track.
  • dEdxError::Float32: error of dEdx.
  • radiusOfInnermostHit::Float32: radius of the innermost hit that has been used in the track fit
  • subdetectorHitNumbers::PVector{Int32}: number of hits in particular subdetectors.Check/set collection variable TrackSubdetectorNames for decoding the indices
  • trackStates::PVector{TrackState}: track states
  • dxQuantities::PVector{Quantity}: different measurements of dx quantities

Relations

  • trackerHits::TrackerHit: hits that have been used to create this track
  • tracks::Track: tracks (segments) that have been combined to create this track

Methods

  • setSubdetectorHitNumbers(object::Track, v::AbstractVector{Int32}): assign a set of values to the subdetectorHitNumbers vector member
  • setTrackStates(object::Track, v::AbstractVector{TrackState}): assign a set of values to the trackStates vector member
  • setDxQuantities(object::Track, v::AbstractVector{Quantity}): assign a set of values to the dxQuantities vector member
  • pushToTrackerHits(obj::Track, robj::TrackerHit): push related object to the trackerHits relation
  • popFromTrackerHits(obj::Track): pop last related object from trackerHits relation
  • pushToTracks(obj::Track, robj::Track): push related object to the tracks relation
  • popFromTracks(obj::Track): pop last related object from tracks relation
EDM4hep.TrackStateType

TrackState

Fields

  • location::Int32: for use with At{Other|IP|FirstHit|LastHit|Calorimeter|Vertex}|LastLocation
  • D0::Float32: transverse impact parameter
  • phi::Float32: azimuthal angle
  • omega::Float32: is the signed curvature of the track in [1/mm].
  • Z0::Float32: longitudinal impact parameter
  • tanLambda::Float32: lambda is the dip angle of the track in r-z
  • time::Float32: time of the track at this trackstate
  • referencePoint::Vector3f: Reference point of the track parameters, e.g. the origin at the IP, or the position of the first/last hits or the entry point into the calorimeter. [mm]
  • covMatrix::SVector{21,Float32}: lower triangular covariance matrix of the track parameters. the order of parameters is d0, phi, omega, z0, tan(lambda), time. the array is a row-major flattening of the matrix.
EDM4hep.TrackerHitType

Tracker hit

  • Author: F.Gaede, DESY

Fields

  • cellID::UInt64: ID of the sensor that created this hit
  • type::Int32: type of raw data hit, either one of edm4hep::RawTimeSeries, edm4hep::SIMTRACKERHIT - see collection parameters "TrackerHitTypeNames" and "TrackerHitTypeValues".
  • quality::Int32: quality bit flag of the hit.
  • time::Float32: time of the hit [ns].
  • eDep::Float32: energy deposited on the hit [GeV].
  • eDepError::Float32: error measured on EDep [GeV].
  • position::Vector3d: hit position in [mm].
  • covMatrix::SVector{6,Float32}: covariance of the position (x,y,z), stored as lower triangle matrix. i.e. cov(x,x) , cov(y,x) , cov(y,y) , cov(z,x) , cov(z,y) , cov(z,z)
  • rawHits::PVector{ObjectID}: raw data hits. Check getType to get actual data type.

Methods

  • setRawHits(object::TrackerHit, v::AbstractVector{ObjectID}): assign a set of values to the rawHits vector member
EDM4hep.TrackerHitPlaneType

Tracker hit plane

  • Author: Placido Fernandez Declara, CERN

Fields

  • cellID::UInt64: ID of the sensor that created this hit
  • type::Int32: type of raw data hit, either one of edm4hep::RawTimeSeries, edm4hep::SIMTRACKERHIT - see collection parameters "TrackerHitTypeNames" and "TrackerHitTypeValues".
  • quality::Int32: quality bit flag of the hit.
  • time::Float32: time of the hit [ns].
  • eDep::Float32: energy deposited on the hit [GeV].
  • eDepError::Float32: error measured on EDep [GeV].
  • u::Vector2f: measurement direction vector, u lies in the x-y plane
  • v::Vector2f: measurement direction vector, v is along z
  • du::Float32: measurement error along the direction
  • dv::Float32: measurement error along the direction
  • position::Vector3d: hit position in [mm].
  • covMatrix::SVector{6,Float32}: covariance of the position (x,y,z), stored as lower triangle matrix. i.e. cov(x,x) , cov(y,x) , cov(y,y) , cov(z,x) , cov(z,y) , cov(z,z)
  • rawHits::PVector{ObjectID}: raw data hits. Check getType to get actual data type.

Methods

  • setRawHits(object::TrackerHitPlane, v::AbstractVector{ObjectID}): assign a set of values to the rawHits vector member
EDM4hep.TrackerPulseType

Reconstructed Tracker Pulse

  • Author: Wenxing Fang, IHEP

Fields

  • cellID::UInt64: cell id.
  • time::Float32: time [ns].
  • charge::Float32: charge [fC].
  • quality::Int16: quality.
  • covMatrix::SVector{3,Float32}: lower triangle covariance matrix of the charge(c) and time(t) measurements.

Relations

  • timeSeries::TimeSeries: Optionally, the timeSeries that has been used to create the pulse can be stored with the pulse.
EDM4hep.Vector4fType

Generic vector for storing classical 4D coordinates in memory. Four momentum helper functions are in edm4hep::utils

Fields

  • x::Float32:
  • y::Float32:
  • z::Float32:
  • t::Float32:
EDM4hep.VertexType

Vertex

  • Author: F.Gaede, DESY

Fields

  • primary::Int32: boolean flag, if vertex is the primary vertex of the event
  • chi2::Float32: chi-squared of the vertex fit
  • probability::Float32: probability of the vertex fit
  • position::Vector3f: [mm] position of the vertex.
  • covMatrix::SVector{6,Float32}: covariance matrix of the position (stored as lower triangle matrix, i.e. cov(xx),cov(y,x),cov(z,x),cov(y,y),... )
  • algorithmType::Int32: type code for the algorithm that has been used to create the vertex - check/set the collection parameters AlgorithmName and AlgorithmType.
  • parameters::PVector{Float32}: additional parameters related to this vertex - check/set the collection parameter "VertexParameterNames" for the parameters meaning.

Relations

  • associatedParticle::POD: reconstructed particle associated to this vertex.

Methods

  • setParameters(object::Vertex, v::AbstractVector{Float32}): assign a set of values to the parameters vector member
EDM4hep.getEDStoreMethod
getEDStore(::Type{ED}, collid::UInt32=0x00000000)

Get the store corresponding to the collid. If it is not specified then obtain a collid frm the data type ED

EDM4hep.hasEDStoreMethod
hasEDStore(collid::UInt32)

Find out if the store with collid is there.

EDM4hep.initEDStoreMethod
initEDStore(::Type{ED}) where ED

Unintialize the store corresponding to type ED

EDM4hep.Histograms.H1DType

H1D(title::String, nbins::Int, min::Float, max::Float, unit::Symbol) Create a 1-dimentional histgram carrying the title and units

EDM4hep.Histograms.H2DType

H2D(title::String, xbins::Int, xmin::Float, xmax::Float, ybins::Int, ymin::Float, ymax::Float, unit::Tuple{Symbol,Symbol}) Create a 2-dimentional histgram carrying the title and units

EDM4hep.Histograms.H3DType

H3D(title::String, xbins::Int, xmin::Float, xmax::Float, ybins::Int, ymin::Float, ymax::Float, zbins::Int, zmin::Float, zmax::Float, unit::Tuple{Symbol,Symbol,Symbol}) Create a 2-dimentional histgram carrying the title and units

EDM4hep.RootIOModule

ROOT I/O module for EDM4hep.jl

It supports both formats: TTree and RNTuple

EDM4hep.RootIO.ReaderType

The Reader structure keeps a reference to the UnROOT LazyTree and caches already built 'layouts' of the EDM4hep types. The layouts maps a set of columns in the LazyTree into an object.

EDM4hep.RootIO.create_getterMethod

create_getter(reader::Reader, bname::String; selection=nothing)

This function creates a getter function for a given branch name. The getter function is a function that takes an event and returns a StructArray with the data of the branch. The getter function is created as a function with the name get_<branchname>. The optional parameter selection is a list of field names to be selected from the branch. If selection is not provided, all fields are selected. The user can use a list of strings, symbols or regular expressions to select the fields.

EDM4hep.RootIO.getMethod

get(reader::Reader, treename::String)

Opens a 'TTree' or 'RNTuple' in the ROOT file (typically the events tree). It returns a 'LazyTree' that allows the user to iterate over events.

EDM4hep.RootIO.getMethod

get(reader::Reader, evt::UnROOT.LazyEvent, bname::String; btype::Type=Any, register=true)

Gets an object collection by its name, with the possibility to overwrite the mapping Julia type or use the type known in the ROOT file (C++ class name). The optonal key parameter register indicates is the collection needs to be registered to the EDStore.