Public Documentation

Documentation for EDM4hep.jl public interface.

Index - Types

Index - Functions

Modules

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

source
EDM4hep.RootIOModule

ROOT I/O module for EDM4hep.jl

It supports both formats: TTree and RNTuple

source
EDM4hep.AnalysisModule

Analysis module for EDM4hep.jl

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

source

Types

This is the list of all types defined for EDM4hep using the PODIO yaml file.

EDM4hep.CaloHitContributionType

Monte Carlo contribution to SimCalorimeterHit

  • Author: EDM4hep authors

Fields

  • PDG::Int32: PDG code of the shower particle that caused this contribution
  • energy::Float32: energy of the this contribution [G]
  • time::Float32: time of this contribution [ns]
  • 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
source
EDM4hep.CalorimeterHitType

Calorimeter hit

  • Author: EDM4hep authors

Fields

  • cellID::UInt64: detector specific (geometrical) cell id
  • energy::Float32: energy of the hit [GeV]
  • energyError::Float32: error of the hit energy [GeV]
  • time::Float32: time of the hit [ns]
  • position::Vector3f: position of the hit in world coordinates [mm]
  • type::Int32: type of hit
source
EDM4hep.ClusterType

Calorimeter Hit Cluster

  • Author: EDM4hep authors

Fields

  • type::Int32: flagword that defines the type of cluster
  • energy::Float32: energy of the cluster [GeV]
  • energyError::Float32: error on the energy [GeV]
  • position::Vector3f: position of the cluster [mm]
  • positionError::CovMatrix3f: covariance matrix of the position
  • iTheta::Float32: Polar angle of the cluster's intrinsic direction (used e.g. for vertexing). Not to be confused with the cluster position seen from IP
  • phi::Float32: Azimuthal angle of the cluster's intrinsic direction (used e.g. for vertexing). Not to be confused with the cluster position seen from IP
  • directionError::Vector3f: covariance matrix of the direction [mm**2]
  • shapeParameters::PVector{Float32}: shape parameters. This should be accompanied by a descriptive list of names in the shapeParameterNames collection level metadata, as a vector of strings with the same ordering
  • subdetectorEnergies::PVector{Float32}: energy observed in a particular subdetector

Relations

  • clusters::Cluster: clusters that have been combined to this cluster
  • hits::CalorimeterHit: hits that have been combined to this cluster

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
source
EDM4hep.CovMatrix2fType

A generic 2 dimensional covariance matrix with values stored in lower triangular form

Fields

  • values::SVector{3,Float32}: the covariance matrix values
source
EDM4hep.CovMatrix3fType

A generic 3 dimensional covariance matrix with values stored in lower triangular form

Fields

  • values::SVector{6,Float32}: the covariance matrix values
source
EDM4hep.CovMatrix4fType

A generic 4 dimensional covariance matrix with values stored in lower triangular form

Fields

  • values::SVector{10,Float32}: the covariance matrix values
source
EDM4hep.CovMatrix6fType

A generic 6 dimensional covariance matrix with values stored in lower triangular form

Fields

  • values::SVector{21,Float32}: the covariance matrix values
source
EDM4hep.EventHeaderType

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

  • Author: EDM4hep authors

Fields

  • eventNumber::UInt64: event number
  • runNumber::UInt32: run number
  • timeStamp::UInt64: time stamp
  • weight::Float64: event weight
  • weights::PVector{Float64}: event weights in case there are multiple. NOTE that weights[0] might not be the same as weight! Event weight names should be stored using the edm4hep::EventWeights name in the file level metadata

Methods

  • setWeights(object::EventHeader, v::AbstractVector{Float64}): assign a set of values to the weights vector member
source
EDM4hep.GeneratorEventParametersType

Generator event parameters

  • Author: EDM4hep authors

Fields

  • eventScale::Float64: event scale
  • alphaQED::Float64: alpha_QED
  • alphaQCD::Float64: alpha_QCD
  • signalProcessId::Int32: id of signal process
  • sqrts::Float64: sqrt(s) [GeV]
  • crossSections::PVector{Float64}: list of cross sections [pb]
  • crossSectionErrors::PVector{Float64}: list of cross section errors [pb]

Relations

  • signalVertex::MCParticle: List of initial state MCParticle that are the source of the hard interaction

Methods

  • setCrossSections(object::GeneratorEventParameters, v::AbstractVector{Float64}): assign a set of values to the crossSections vector member
  • setCrossSectionErrors(object::GeneratorEventParameters, v::AbstractVector{Float64}): assign a set of values to the crossSectionErrors vector member
  • pushToSignalVertex(obj::GeneratorEventParameters, robj::MCParticle): push related object to the signalVertex relation
  • popFromSignalVertex(obj::GeneratorEventParameters): pop last related object from signalVertex relation
source
EDM4hep.GeneratorPdfInfoType

Generator pdf information

  • Author: EDM4hep authors

Fields

  • partonId::SVector{2,Int32}: Parton PDG id
  • lhapdfId::SVector{2,Int32}: LHAPDF PDF id (see https://lhapdf.hepforge.org/pdfsets.html)
  • x::SVector{2,Float64}: Parton momentum fraction
  • xf::SVector{2,Float64}: PDF value
  • scale::Float64: Factorisation scale [GeV]
source
EDM4hep.MCParticleType

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

  • Author: EDM4hep authors

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 wrt. the event, e.g. for preassigned decays or decays in flight from the simulator [ns]
  • mass::Float64: mass of the particle [GeV]
  • vertex::Vector3d: production vertex of the particle [mm]
  • endpoint::Vector3d: endpoint of the particle [mm]
  • momentum::Vector3d: particle 3-momentum at the production vertex [GeV]
  • momentumAtEndpoint::Vector3d: particle 3-momentum at the endpoint [GeV]
  • spin::Vector3f: spin (helicity) vector of the particle

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
source
EDM4hep.ParticleIDType

ParticleID

  • Author: EDM4hep authors

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

Relations

  • particle::ReconstructedParticle: the particle from which this PID has been computed

Methods

  • setParameters(object::ParticleID, v::AbstractVector{Float32}): assign a set of values to the parameters vector member
source
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
source
EDM4hep.RawCalorimeterHitType

Raw calorimeter hit

  • Author: EDM4hep authors

Fields

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

Raw data of a detector readout

  • Author: EDM4hep authors

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
source
EDM4hep.RecDqdxType

dN/dx or dE/dx info of a Track

  • Author: EDM4hep authors

Fields

  • dQdx::Quantity: the reconstructed dEdx or dNdx and its error

Relations

  • track::Track: the corresponding track
source
EDM4hep.ReconstructedParticleType

Reconstructed Particle

  • Author: EDM4hep authors

Fields

  • PDG::Int32: PDG of the reconstructed particle.
  • energy::Float32: energy of the reconstructed particle. Four momentum state is not kept consistent internally [GeV]
  • momentum::Vector3f: particle momentum. Four momentum state is not kept consistent internally [GeV]
  • referencePoint::Vector3f: reference, i.e. where the particle has been measured [mm]
  • charge::Float32: charge of the reconstructed particle
  • mass::Float32: mass of the reconstructed particle, set independently from four vector. Four momentum state is not kept consistent internally [GeV]
  • goodnessOfPID::Float32: overall goodness of the PID on a scale of [0;1]
  • covMatrix::CovMatrix4f: covariance matrix of the reconstructed particle 4vector

Relations

  • decayVertex::Vertex: decay vertex for the particle (if it is a composite 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

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
source
EDM4hep.SimCalorimeterHitType

Simulated calorimeter hit

  • Author: EDM4hep authors

Fields

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

Relations

  • contributions::CaloHitContribution: Monte Carlo step contributions

Methods

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

Simulated tracker hit

  • Author: EDM4hep authors

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 [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 [mm]
  • momentum::Vector3f: the 3-momentum of the particle at the hits position [GeV]

Relations

  • particle::MCParticle: MCParticle that caused the hit
source
EDM4hep.TimeSeriesType

Calibrated Detector Data

  • Author: EDM4hep authors

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
source
EDM4hep.TrackType

Reconstructed track

  • Author: EDM4hep authors

Fields

  • type::Int32: flagword that defines the type of track
  • chi2::Float32: Chi^2 of the track fit
  • ndf::Int32: number of degrees of freedom of the track fit
  • Nholes::Int32: number of holes on track
  • subdetectorHitNumbers::PVector{Int32}: number of hits in particular subdetectors
  • subdetectorHoleNumbers::PVector{Int32}: number of holes in particular subdetectors
  • trackStates::PVector{TrackState}: track states

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
  • setSubdetectorHoleNumbers(object::Track, v::AbstractVector{Int32}): assign a set of values to the subdetectorHoleNumbers vector member
  • setTrackStates(object::Track, v::AbstractVector{TrackState}): assign a set of values to the trackStates 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
source
EDM4hep.TrackStateType

Parametrized description of a particle track

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 [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 [ns]
  • 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::CovMatrix6f: covariance matrix of the track parameters.
source
EDM4hep.TrackerHitType

Tracker hit interface class

  • Author: Thomas Madlener, DESY

Fields

  • cellID::UInt64: ID of the sensor that created this hit
  • type::Int32: type of the raw data hit
  • 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 [mm]
source
EDM4hep.TrackerHit3DType

Tracker hit

  • Author: EDM4hep authors

Fields

  • cellID::UInt64: ID of the sensor that created this hit
  • type::Int32: type of raw data hit
  • 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 [mm]
  • covMatrix::CovMatrix3f: covariance matrix of the position (x,y,z)
source
EDM4hep.TrackerHitPlaneType

Tracker hit plane

  • Author: EDM4hep authors

Fields

  • cellID::UInt64: ID of the sensor that created this hit
  • type::Int32: type of raw data hit
  • 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 [mm]
  • covMatrix::CovMatrix3f: covariance of the position (x,y,z)
source
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:
source
EDM4hep.VertexType

Vertex

  • Author: EDM4hep authors

Fields

  • type::UInt32: Flagword that defines the type of the vertex, see reserved bits for more information
  • chi2::Float32: chi-squared of the vertex fit
  • ndf::Int32: number of degrees of freedom of the vertex fit
  • position::Vector3f: [mm] position of the vertex
  • covMatrix::CovMatrix3f: covariance matrix of the position
  • algorithmType::Int32: type code for the algorithm that has been used to create the vertex
  • parameters::PVector{Float32}: additional parameters related to this vertex

Relations

  • particles::POD: particles that have been used to form this vertex, aka the decay particles emerging from this vertex

Methods

  • setParameters(object::Vertex, v::AbstractVector{Float32}): assign a set of values to the parameters vector member
  • pushToParticles(obj::Vertex, robj::POD): push related object to the particles relation
  • popFromParticles(obj::Vertex): pop last related object from particles relation
source
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.

source
EDM4hep.Histograms.H1DType

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

source
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-dimensional histogram carrying the title and units.

source
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-dimensional histogram carrying the title and units.

source

Functions

EDM4hep.getEDCollectionMethod
getEDCollection(::Type{ED}, collid::UInt32=0x00000000)

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

source
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.

source
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.

source
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 optional key parameter register indicates if the collection needs to be registered to the EDStore.

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

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

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

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

source
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 different tasks running on different threads if the option argument mt is set to true. The results in data for all the chunks will be merged at the end of the analysis.

source