ATLAS Offline Software
Loading...
Searching...
No Matches
egammaMVAFunctions.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef EGAMMAMVACALIB_EGAMMAMVAFUNCTIONS
6#define EGAMMAMVACALIB_EGAMMAMVAFUNCTIONS
7
8#include "xAODEgamma/Egamma.h"
10#include "xAODEgamma/Photon.h"
14#include "egammaMVALayerDepth.h"
16
17#include "TLorentzVector.h"
18
19#include <functional>
20#include <string>
21#include <unordered_map>
22#include <cmath>
23#include <memory>
24#include <stdexcept>
25
26// for the ConversionHelper (deprecated?)
28
39
40// Changing the definition of the functions means breaking backward
41// compatibility with previous version of the MVA calibrations.
42
44{
45 // inline functions to avoid duplicates problem during linking (and performance)
46 // cluster functions
47 // REMEMBER to add the functions using corrected layer energies
48 inline float compute_cl_eta(const xAOD::CaloCluster& cluster) { return cluster.eta(); }
49 inline float compute_cl_phi(const xAOD::CaloCluster& cluster) { return cluster.phi(); }
50 inline float compute_cl_e(const xAOD::CaloCluster& cluster) { return cluster.e(); }
51 inline float compute_cl_etaCalo(const xAOD::CaloCluster& cluster) {
52 double tmp = 0.;
54 throw std::runtime_error("etaCalo not found in CaloCluster object");
55 }
56 return tmp;
57 }
58 inline float compute_cl_phiCalo(const xAOD::CaloCluster& cluster) {
59 double tmp = 0.;
61 throw std::runtime_error("phiCalo not found in CaloCluster object");
62 }
63 return tmp;
64 }
65 inline float compute_cl_etas1(const xAOD::CaloCluster& cluster) { return cluster.etaBE(1); }
66 inline float compute_cl_etas2(const xAOD::CaloCluster& cluster) { return cluster.etaBE(2); }
67 inline float compute_rawcl_Es0(const xAOD::CaloCluster& cl) { return cl.energyBE(0); }
68 /*inline std::function<float(const xAOD::CaloCluster&)> compute_rawcl_Es0_auto(bool use_corrected)
69 {
70 if (use_corrected) return [](const xAOD::CaloCluster& cl) { return cl.energyBE(0); };
71 else return [](const xAOD::CaloCluster& cl) { return cl.energyBE(0); };
72 }*/
73 inline float compute_rawcl_Es1(const xAOD::CaloCluster& cl) { return cl.energyBE(1); }
74 inline float compute_rawcl_Es2(const xAOD::CaloCluster& cl) { return cl.energyBE(2); }
75 inline float compute_rawcl_Es3(const xAOD::CaloCluster& cl) { return cl.energyBE(3); }
76
78 static const SG::ConstAccessor<double> acc ("correctedcl_Es0");
79 return acc.isAvailable(cl) ? acc(cl) : cl.energyBE(0);
80 }
82 static const SG::ConstAccessor<double> acc ("correctedcl_Es1");
83 return acc.isAvailable(cl) ? acc(cl) : cl.energyBE(1);
84 }
86 static const SG::ConstAccessor<double> acc ("correctedcl_Es2");
87 return acc.isAvailable(cl) ? acc(cl) : cl.energyBE(2);
88 }
90 static const SG::ConstAccessor<double> acc ("correctedcl_Es3");
91 return acc.isAvailable(cl) ? acc(cl) : cl.energyBE(3);
92 }
93
94 inline float compute_rawcl_Eacc(const xAOD::CaloCluster& cl) { return cl.energyBE(1) + cl.energyBE(2) + cl.energyBE(3); }
95 inline float compute_rawcl_f0(const xAOD::CaloCluster& cl) { return cl.energyBE(0) / (cl.energyBE(1) + cl.energyBE(2) + cl.energyBE(3)); }
96
99
100
101 inline float compute_calibHitsShowerDepth(const std::array<float, 4>& cl, float eta)
102 {
103 const float denominator = cl[0] + cl[1] + cl[2] + cl[3];
104 if (denominator == 0) return 0.;
105
106 const std::array<float, 4> radius(get_MVAradius(eta));
107
108 // loop unrolling
109 return (radius[0] * cl[0]
110 + radius[1] * cl[1]
111 + radius[2] * cl[2]
112 + radius[3] * cl[3]) / denominator;
113 }
114
116 {
117 const std::array<float, 4> cluster_array {{ compute_rawcl_Es0(cl),
120 compute_rawcl_Es3(cl) }};
121 return compute_calibHitsShowerDepth(cluster_array, compute_cl_eta(cl));
122 }
123
125 const std::array<float, 4> cluster_array {{ compute_correctedcl_Es0(cl),
129 return compute_calibHitsShowerDepth(cluster_array, compute_cl_eta(cl));
130 }
131
132 // ------------------------------------------------------------------
133 // Forward-electron getters
134 // ------------------------------------------------------------------
135
136 inline float compute_et(const xAOD::CaloCluster& cl) {
137 const float e = cl.e();
138 const float eta = cl.eta();
139 const float c = std::cosh(eta);
140 return (c != 0.f) ? (e / c) : 0.f;
141 }
142 // be 100% sure what variables the FE BDT should use. For the run 2
143 inline float cl_getMoment(const xAOD::CaloCluster& cl,
145 const char* name) {
146 double tmp = 0.;
147 if (!cl.retrieveMoment(m, tmp)) {
148 throw std::runtime_error(std::string("Forward-electron missing moment: ") + name);
149 }
150 return static_cast<float>(tmp);
151 }
152
153 inline float compute_cl_significance (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::SIGNIFICANCE , "SIGNIFICANCE"); }
154 inline float compute_cl_secondLambda (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::SECOND_LAMBDA , "SECOND_LAMBDA"); }
155 inline float compute_cl_lateral (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::LATERAL , "LATERAL"); }
156 inline float compute_cl_longitudinal (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::LONGITUDINAL , "LONGITUDINAL"); }
157 inline float compute_cl_fracMax (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::ENG_FRAC_MAX , "ENG_FRAC_MAX"); }
158
159 inline float compute_cl_secondR (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::SECOND_R , "SECOND_R"); }
160 inline float compute_cl_centerLambda (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::CENTER_LAMBDA , "CENTER_LAMBDA"); }
161 inline float compute_cl_secondDensity(const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::SECOND_ENG_DENS,"SECOND_ENG_DENS"); }
162 inline float compute_cl_x (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::CENTER_X , "CENTER_X"); }
163 inline float compute_cl_y (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::CENTER_Y , "CENTER_Y"); }
164 inline float compute_cl_z (const xAOD::CaloCluster& cl){ return cl_getMoment(cl, xAOD::CaloCluster::CENTER_Z , "CENTER_Z"); }
165
166 inline float compute_cl_secondR_fudge(const xAOD::Egamma& eg) {
167 static const SG::AuxElement::Accessor<float> accR2("SECOND_R");
168 if (accR2.isAvailable(eg)) { return accR2(eg); }
169 return -1.;
170 }
171
172 inline float compute_eta_FCAL(const xAOD::CaloCluster& cl){
173 float x = compute_cl_x(cl);
174 float y = compute_cl_y(cl);
175 float z = compute_cl_z(cl);
176 float theta = std::acos(z/std::sqrt(x*x+y*y+z*z));
177 return -std::log(std::tan(theta/2.));
178 }
179 inline float compute_etaMod_FCAL(const xAOD::CaloCluster& cl){
180 return std::fmod(std::abs(compute_eta_FCAL(cl)),0.15);
181 }
183 return std::floor(std::abs(compute_eta_FCAL(cl))/0.15);
184 }
185 inline float compute_etaMod_EMEC(const xAOD::CaloCluster& cl){
186 return std::fmod(std::abs(cl.eta()),0.1);
187 }
189 return std::floor(std::abs(cl.eta())/0.1);
190 }
191 inline float compute_phiMod_EMEC(const xAOD::CaloCluster& cl){
192 static const float cz = std::numbers::pi/16.;
193 float phi_mod = std::fmod(cl.phi(), cz);
194 if (phi_mod < 0) phi_mod += cz;
195 return phi_mod;
196 }
197 inline float compute_R12_EMEC(const xAOD::CaloCluster& cl){
198 double es1 = cl.energy_max(CaloSampling::EME2);
199 double es2 = cl.energy_max(CaloSampling::EME3);
200 return es2 != 0 ? float(es1/es2) : -1;
201 }
202
203 // electron functions
204 inline float compute_el_charge(const xAOD::Electron& el) { return el.charge(); }
205 inline float compute_el_tracketa(const xAOD::Electron& el) { return el.trackParticle()->eta(); }
206 inline float compute_el_trackpt(const xAOD::Electron& el) { return el.trackParticle()->pt(); }
207 inline float compute_el_trackz0(const xAOD::Electron& el) { return el.trackParticle()->z0(); }
208 inline float compute_el_refittedTrack_qoverp(const xAOD::Electron& el) { return el.trackParticle()->qOverP(); }
209 inline int compute_el_author(const xAOD::Electron& el) {
210 static const SG::ConstAccessor<unsigned short int> acc ("author");
211 return acc (el);
212 }
213
214 // photon functions
215 inline int compute_ph_convFlag(const xAOD::Photon& ph) {
216 const auto original = xAOD::EgammaHelpers::conversionType(&ph);
217 if (original == 3) return 2;
218 else if (original != 0) return 1;
219 else return original;
220 }
221
222 // a utility function
224 {
225 if (!tp) return 0;
226 for (unsigned int i = 0; i < tp->numberOfParameters(); ++i) {
228 return hypot(tp->parameterPX(i), tp->parameterPY(i));
229 }
230 }
231 return tp->pt();
232 }
233
234 // define a few without using conversion helper
235
237 inline float compute_ptconv_decor(const xAOD::Photon* ph)
238 {
239 static const SG::AuxElement::Accessor<float> accPx("px");
240 static const SG::AuxElement::Accessor<float> accPy("py");
241
242 auto vx = ph->vertex();
243 return vx ? std::hypot(accPx(*vx), accPy(*vx)) : 0.0;
244 }
245
247 inline float compute_ptconv(const xAOD::Photon* ph)
248 {
249 auto vx = ph->vertex();
250 if (!vx) return 0.0;
251
252 TLorentzVector sum;
253 if (vx->trackParticle(0)) sum += vx->trackParticle(0)->p4();
254 if (vx->trackParticle(1)) sum += vx->trackParticle(1)->p4();
255 return sum.Perp();
256 }
257
258 inline float compute_pt1conv(const xAOD::Photon* ph)
259 {
260 static const SG::AuxElement::Accessor<float> accPt1("pt1");
261
262 const xAOD::Vertex* vx = ph->vertex();
263 if (!vx) return 0.0;
264 if (accPt1.isAvailable(*vx)) {
265 return accPt1(*vx);
266 } else {
268 }
269 }
270
271 inline float compute_pt2conv(const xAOD::Photon* ph)
272 {
273 static const SG::AuxElement::Accessor<float> accPt2("pt2");
274
275 const xAOD::Vertex* vx = ph->vertex();
276 if (!vx) return 0.0;
277 if (accPt2.isAvailable(*vx)) {
278 return accPt2(*vx);
279 } else {
281 }
282 }
283
284 // using template to avoid rewriting code for 1st, 2nd track and
285 // for all the summary types
286 template<int itrack, xAOD::SummaryType summary>
287 inline int compute_convtrkXhits(const xAOD::Photon* ph) {
288 const auto vx = ph->vertex();
289 if (!vx) return 0.;
290
291 if (vx->trackParticle(0)) {
292 uint8_t hits;
293 if (vx->trackParticle(itrack)->summaryValue(hits, summary)) {
294 return hits;
295 }
296 }
297 return 0.;
298 }
299
304
305 // The functions to return the dictionaries of functions,
306 // i.e., the variable name to function
307
309 using funcMap_t = std::unordered_map<std::string,
310 std::function<float(const xAOD::Egamma*, const xAOD::CaloCluster*)> >;
311
313 std::unique_ptr<funcMap_t> initializeElectronFuncs(bool useLayerCorrected);
314
316 std::unique_ptr<funcMap_t> initializeUnconvertedPhotonFuncs(bool useLayerCorrected);
317
319 std::unique_ptr<funcMap_t> initializeConvertedPhotonFuncs(bool useLayerCorrected);
320
322 std::unique_ptr<funcMap_t> initializeForwardElectronFuncs(bool useLayerCorrected);
323
327 {
329 : asg::AsgMessaging("ConversionHelper"),
330 m_vertex(ph ? ph->vertex() : nullptr),
331 m_tp0(m_vertex ? m_vertex->trackParticle(0) : nullptr),
332 m_tp1(m_vertex ? m_vertex->trackParticle(1) : nullptr),
333 m_pt1conv(0.), m_pt2conv(0.)
334 {
335
336 ATH_MSG_DEBUG("init conversion helper");
337 if (!m_vertex) return;
338
339 static const SG::AuxElement::Accessor<float> accPt1("pt1");
340 static const SG::AuxElement::Accessor<float> accPt2("pt2");
341 if (accPt1.isAvailable(*m_vertex) && accPt2.isAvailable(*m_vertex))
342 {
343 m_pt1conv = accPt1(*m_vertex);
344 m_pt2conv = accPt2(*m_vertex);
345 }
346 else
347 {
348 ATH_MSG_WARNING("pt1/pt2 not available, will approximate from first measurements");
351 }
352 }
353
355 : ConversionHelper(&ph) { } // delegating constr
356
357 inline float ph_Rconv() const { return m_vertex ? hypot(m_vertex->position().x(), m_vertex->position().y()) : 0; }
358 inline float ph_zconv() const { return m_vertex ? m_vertex->position().z() : 0.; }
359 inline int ph_convtrk1nPixHits() const {
360 if (!m_tp0) { return 0; }
361 uint8_t hits = 0;
362 if (m_tp0->summaryValue(hits, xAOD::numberOfPixelHits)) { return hits; }
363 else {
364 ATH_MSG_WARNING("cannot read xAOD::numberOfPixelHits");
365 return 0;
366 }
367 }
368 inline int ph_convtrk2nPixHits() const {
369 if (!m_tp1) return 0;
370 uint8_t hits;
371 if (m_tp1->summaryValue(hits, xAOD::numberOfPixelHits)) { return hits; }
372 else {
373 ATH_MSG_WARNING("cannot read xAOD::numberOfPixelHits");
374 return 0;
375 }
376 }
377 inline int ph_convtrk1nSCTHits() const {
378 if (!m_tp0) { return 0; }
379 uint8_t hits;
380 if (m_tp0->summaryValue(hits, xAOD::numberOfSCTHits)) { return hits; }
381 else {
382 ATH_MSG_WARNING("cannot read xAOD::numberOfSCTHits");
383 return 0;
384 }
385 }
386 inline int ph_convtrk2nSCTHits() const {
387 if (!m_tp1) { return 0; }
388 uint8_t hits;
389 if (m_tp1->summaryValue(hits, xAOD::numberOfSCTHits)) { return hits; }
390 else {
391 ATH_MSG_WARNING("cannot read xAOD::numberOfSCTHits");
392 return 0;
393 }
394 }
395 inline float ph_pt1conv() const { return m_pt1conv; }
396 inline float ph_pt2conv() const { return m_pt2conv; }
397 inline float ph_ptconv() const {
398 // TODO: evaluate if move to this new definition, by now keep the previous one
399 // to be consistent with the training
400 // return m_vertex ? xAOD::EgammaHelpers::momentumAtVertex(*m_vertex).perp() : 0.;
401 TLorentzVector sum;
402 if (m_tp0) sum += m_tp0->p4();
403 if (m_tp1) sum += m_tp1->p4();
404 return sum.Perp();
405 }
406 private:
411 };
412
413
414} // end namespace
415
416#endif
Scalar eta() const
pseudorapidity method
Scalar theta() const
theta method
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Helper class to provide constant type-safe access to aux data.
#define y
#define x
#define z
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:573
Helper class to provide constant type-safe access to aux data.
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
Class mimicking the AthMessaging class from the offline software.
AsgMessaging(const std::string &name)
Constructor with a name.
bool retrieveMoment(MomentType type, double &value) const
Retrieve individual moment.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double e() const
The total energy of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
MomentType
Enums to identify different moments.
@ SECOND_ENG_DENS
Second Moment in E/V.
@ SECOND_LAMBDA
Second Moment in .
@ ETACALOFRAME
Eta in the calo frame (for egamma)
@ LATERAL
Normalized lateral moment.
@ LONGITUDINAL
Normalized longitudinal moment.
@ ENG_FRAC_MAX
Energy fraction of hottest cell.
@ SECOND_R
Second Moment in .
@ PHICALOFRAME
Phi in the calo frame (for egamma)
@ CENTER_LAMBDA
Shower depth at Cluster Centroid.
@ SIGNIFICANCE
Cluster significance.
@ CENTER_Z
Cluster Centroid ( )
@ CENTER_X
Cluster Centroid ( )
@ CENTER_Y
Cluster Centroid ( )
float etaBE(const unsigned layer) const
Get the eta in one layer of the EM Calo.
const xAOD::Vertex * vertex(size_t index=0) const
Pointer to the xAOD::Vertex/es that match the photon candidate.
Definition Photon_v1.cxx:61
float parameterPX(unsigned int index) const
Returns the parameter x momentum component, for 'index'.
float parameterPY(unsigned int index) const
Returns the parameter y momentum component, for 'index'.
xAOD::ParameterPosition parameterPosition(unsigned int index) const
Return the ParameterPosition of the parameters at 'index'.
size_t numberOfParameters() const
Returns the number of additional parameters stored in the TrackParticle.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
const TrackParticle * trackParticle(size_t i) const
Get the pointer to a given track that was used in vertex reco.
std::array< float, 4 > get_MVAradius(float eta)
helper function to compute shower depth
These functions are for calculating variables used by the MVA calibration.
float compute_rawcl_calibHitsShowerDepth(const xAOD::CaloCluster &cl)
float compute_correctedcl_f0(const xAOD::CaloCluster &cl)
float compute_cl_phiCalo(const xAOD::CaloCluster &cluster)
float compute_ptconv(const xAOD::Photon *ph)
This ptconv is the old one used by MVACalib.
std::unique_ptr< funcMap_t > initializeUnconvertedPhotonFuncs(bool useLayerCorrected)
A function to build the map for uncoverted photons.
float compute_cl_y(const xAOD::CaloCluster &cl)
float compute_correctedcl_Eacc(const xAOD::CaloCluster &cl)
float compute_cl_longitudinal(const xAOD::CaloCluster &cl)
float compute_cl_phi(const xAOD::CaloCluster &cluster)
float compute_el_charge(const xAOD::Electron &el)
float compute_rawcl_Es1(const xAOD::CaloCluster &cl)
float compute_etaMod_FCAL(const xAOD::CaloCluster &cl)
float compute_cl_z(const xAOD::CaloCluster &cl)
std::unique_ptr< funcMap_t > initializeElectronFuncs(bool useLayerCorrected)
A function to build the map for electrons.
float compute_correctedcl_Es3(const xAOD::CaloCluster &cl)
float compute_cellIndex_EMEC(const xAOD::CaloCluster &cl)
float compute_cl_secondR(const xAOD::CaloCluster &cl)
float compute_cl_centerLambda(const xAOD::CaloCluster &cl)
float compute_rawcl_Eacc(const xAOD::CaloCluster &cl)
int compute_convtrkXhits(const xAOD::Photon *ph)
float compute_cl_etas2(const xAOD::CaloCluster &cluster)
float compute_el_trackpt(const xAOD::Electron &el)
std::unique_ptr< funcMap_t > initializeForwardElectronFuncs(bool useLayerCorrected)
NEW: A function to build the map for forward electrons.
float compute_cl_etaCalo(const xAOD::CaloCluster &cluster)
float compute_cl_e(const xAOD::CaloCluster &cluster)
float compute_rawcl_Es2(const xAOD::CaloCluster &cl)
float compute_correctedcl_Es1(const xAOD::CaloCluster &cl)
float compute_pt2conv(const xAOD::Photon *ph)
float compute_rawcl_Es0(const xAOD::CaloCluster &cl)
float compute_cl_secondR_fudge(const xAOD::Egamma &eg)
float compute_cl_lateral(const xAOD::CaloCluster &cl)
float compute_calibHitsShowerDepth(const std::array< float, 4 > &cl, float eta)
float compute_rawcl_Es3(const xAOD::CaloCluster &cl)
float compute_ptconv_decor(const xAOD::Photon *ph)
This ptconv function uses the vertex decorations.
float cl_getMoment(const xAOD::CaloCluster &cl, xAOD::CaloCluster::MomentType m, const char *name)
float compute_el_tracketa(const xAOD::Electron &el)
int compute_el_author(const xAOD::Electron &el)
float compute_cellIndex_FCAL(const xAOD::CaloCluster &cl)
float compute_correctedcl_Es2(const xAOD::CaloCluster &cl)
float compute_et(const xAOD::CaloCluster &cl)
float compute_cl_eta(const xAOD::CaloCluster &cluster)
int compute_convtrk2nSCTHits(const xAOD::Photon *ph)
int compute_convtrk1nPixHits(const xAOD::Photon *ph)
float compute_cl_secondLambda(const xAOD::CaloCluster &cl)
float compute_el_refittedTrack_qoverp(const xAOD::Electron &el)
float compute_R12_EMEC(const xAOD::CaloCluster &cl)
float compute_etaMod_EMEC(const xAOD::CaloCluster &cl)
int compute_ph_convFlag(const xAOD::Photon &ph)
std::unique_ptr< funcMap_t > initializeConvertedPhotonFuncs(bool useLayerCorrected)
A function to build the map for converted photons.
float compute_eta_FCAL(const xAOD::CaloCluster &cl)
float compute_correctedcl_Es0(const xAOD::CaloCluster &cl)
float getPtAtFirstMeasurement(const xAOD::TrackParticle *tp)
float compute_el_trackz0(const xAOD::Electron &el)
float compute_cl_etas1(const xAOD::CaloCluster &cluster)
float compute_cl_x(const xAOD::CaloCluster &cl)
float compute_cl_significance(const xAOD::CaloCluster &cl)
float compute_cl_secondDensity(const xAOD::CaloCluster &cl)
float compute_correctedcl_calibHitsShowerDepth(const xAOD::CaloCluster &cl)
float compute_phiMod_EMEC(const xAOD::CaloCluster &cl)
float compute_rawcl_f0(const xAOD::CaloCluster &cl)
std::unordered_map< std::string, std::function< float(const xAOD::Egamma *, const xAOD::CaloCluster *)> > funcMap_t
Define the map type since it's long.
int compute_convtrk2nPixHits(const xAOD::Photon *ph)
float compute_cl_fracMax(const xAOD::CaloCluster &cl)
float compute_pt1conv(const xAOD::Photon *ph)
int compute_convtrk1nSCTHits(const xAOD::Photon *ph)
xAOD::EgammaParameters::ConversionType conversionType(const xAOD::Photon *ph)
return the photon conversion type (see EgammaEnums)
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.
Egamma_v1 Egamma
Definition of the current "egamma version".
Definition Egamma.h:17
Photon_v1 Photon
Definition of the current "egamma version".
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
Electron_v1 Electron
Definition of the current "egamma version".
@ FirstMeasurement
Parameter defined at the position of the 1st measurement.