ATLAS Offline Software
Loading...
Searching...
No Matches
DerivationFramework::EGElectronAmbiguityTool Class Reference

#include <EGElectronAmbiguityTool.h>

Inheritance diagram for DerivationFramework::EGElectronAmbiguityTool:
Collaboration diagram for DerivationFramework::EGElectronAmbiguityTool:

Classes

struct  DecorHandles

Public Member Functions

virtual StatusCode initialize () override final
virtual StatusCode addBranches (const EventContext &ctx) const override final

Private Member Functions

StatusCode decorateSimple (DecorHandles &dh, std::unique_ptr< ConstDataVector< xAOD::TrackParticleContainer > > &tpC, const xAOD::Electron *ele, const xAOD::Vertex *pvtx) const

Private Attributes

SG::ReadHandleKey< xAOD::ElectronContainerm_containerName
SG::ReadHandleKey< xAOD::VertexContainerm_VtxContainerName
SG::ReadHandleKey< xAOD::TrackParticleContainerm_tpContainerName
SG::ReadHandleKey< xAOD::TrackParticleContainerm_tpCName
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_drv { this, "DFCommonSimpleConvRadius", m_containerName, "DFCommonSimpleConvRadius", "" }
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_dphiv { this, "DFCommonSimpleConvPhi", m_containerName, "DFCommonSimpleConvPhi", "" }
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_dmee { this, "DFCommonSimpleMee", m_containerName, "DFCommonSimpleMee", "" }
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_dmeeVtx { this, "DFCommonSimpleMeeAtVtx", m_containerName, "DFCommonSimpleMeeAtVtx", "" }
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_dsep { this, "DFCommonSimpleSeparation", m_containerName, "DFCommonSimpleSeparation", "" }
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_dambi { this, "DFCommonAddAmbiguity", m_containerName, "DFCommonAddAmbiguity", "" }
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_dtrv { this, "DFCommonProdTrueRadius", m_containerName, "DFCommonProdTrueRadius", "" }
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_dtpv { this, "DFCommonProdTruePhi", m_containerName, "DFCommonProdTruePhi", "" }
SG::WriteDecorHandleKey< xAOD::ElectronContainerm_dtzv { this, "DFCommonProdTrueZ", m_containerName, "DFCommonProdTrueZ", "" }
Gaudi::Property< bool > m_isMC {this, "isMC", false}
Gaudi::Property< double > m_elepTCut
Gaudi::Property< std::string > m_idCut
Gaudi::Property< unsigned int > m_nSiCut
Gaudi::Property< double > m_dctCut
Gaudi::Property< double > m_sepCut
Gaudi::Property< double > m_dzCut
Gaudi::Property< double > m_rvECCut
Gaudi::Property< double > m_meeAtVtxECCut
Gaudi::Property< double > m_meeICCut

Friends

struct DecorHandles

Detailed Description

Definition at line 29 of file EGElectronAmbiguityTool.h.

Member Function Documentation

◆ addBranches()

StatusCode DerivationFramework::EGElectronAmbiguityTool::addBranches ( const EventContext & ctx) const
finaloverridevirtual

Definition at line 83 of file EGElectronAmbiguityTool.cxx.

84{
85
86 DecorHandles dh (*this, ctx);
87
88 static const SG::AuxElement::ConstAccessor<char> aidCut(m_idCut);
89
90 // retrieve primary vertex
91 const xAOD::Vertex* pvtx(nullptr);
92 SG::ReadHandle<xAOD::VertexContainer> inputVtx{ m_VtxContainerName, ctx };
93 if (inputVtx.ptr()){
94 const xAOD::VertexContainer* vtxC = inputVtx.ptr();
95 for (const auto* vertex : *vtxC) {
96 if (vertex->vertexType() == xAOD::VxType::VertexType::PriVtx) {
97 pvtx = vertex;
98 break;
99 }
100 }
101 }
102
103 // retrieve electron container
104 SG::ReadHandle<xAOD::ElectronContainer> inputElectrons{ m_containerName,
105 ctx };
106 const xAOD::ElectronContainer* eleC = inputElectrons.ptr();
107
108 if (!eleC) {
110 "Couldn't retrieve Electron container with key: " << m_containerName);
111 return StatusCode::FAILURE;
112 }
113 if (!pvtx) {
114 ATH_MSG_DEBUG("No primary vertex found. Setting default values.");
115 for (const xAOD::Electron* iele : *eleC) {
116 dh.drv(*iele) = -1;
117 dh.dphiv(*iele) = -1;
118 dh.dmee(*iele) = -1;
119 dh.dmeeVtx(*iele) = -1;
120 dh.dsep(*iele) = -1;
121 dh.dambi(*iele) = -1;
122 dh.dtrv(*iele) = -1;
123 dh.dtpv(*iele) = -1;
124 dh.dtzv(*iele) = -1;
125 }
126 return StatusCode::SUCCESS;
127 }
128 ATH_MSG_DEBUG("Pvx z = " << pvtx->z() << ", number of electrons "
129 << eleC->size());
130
131 // Make a container of selected tracks : with Si hits, close to electron track
132 SG::ReadHandle<xAOD::TrackParticleContainer> inputTrackParticles{
134 };
135 const xAOD::TrackParticleContainer* idtpC = inputTrackParticles.ptr();
136
137 std::set<const xAOD::TrackParticle*> alreadyStored;
138 std::set<const xAOD::TrackParticle*> eleIDtpStored, eleGSFtpStored;
139 auto closeByTracks =
140 std::make_unique<ConstDataVector<xAOD::TrackParticleContainer>>(
142
143 for (const auto* ele : *eleC) {
144
145 dh.dambi(*ele) = -1;
146
147 // Electron preselection
148 if (ele->pt() < m_elepTCut ||
149 m_idCut.empty() || !aidCut.isAvailable(*ele) || !aidCut(*ele))
150 continue;
151
152 // Just for debug
153 const xAOD::TrackParticle* eleGSFtp = ele->trackParticle();
154 eleGSFtpStored.insert(eleGSFtp);
155
156 const xAOD::TrackParticle* eleIDtp =
158 eleIDtpStored.insert(eleIDtp);
159
160 // The loop on track
161 for (const auto* tp : *idtpC) {
162
163 // Keep the electron track (I build a container to run vertexing on it...)
164 if (tp == eleIDtp) {
165 closeByTracks->push_back(tp);
166 alreadyStored.insert(tp);
167 continue;
168 }
169
170 // potential candidate to store if not already there
171 if (alreadyStored.find(tp) != alreadyStored.end())
172 continue;
173
175 // if (tp->charge() * ele->charge() > 0)
176 // continue;
177
178 // Close-by
179 double dR = eleIDtp->p4().DeltaR(tp->p4());
180 double dz = std::abs(eleIDtp->z0() - tp->z0()) * sin(eleIDtp->theta());
181 if (dR >= 0.3 || dz >= m_dzCut)
182 continue;
183
184 // With minimum number of Si hits
186 continue;
187
188 alreadyStored.insert(tp);
189
190 closeByTracks->push_back(tp);
191 }
192 }
193
194 if (closeByTracks->empty())
195 return StatusCode::SUCCESS;
196
197 if (msgLvl(MSG::DEBUG)) {
198 SG::ReadHandle<xAOD::TrackParticleContainer> tpCReadHandle{ m_tpCName,
199 ctx };
200 const xAOD::TrackParticleContainer* tpC = tpCReadHandle.ptr();
201
202 ATH_MSG_DEBUG("Number of input tracks "
203 << idtpC->size() << " , number of selected close-by tracks "
204 << closeByTracks->size() << " , number of GSF tracks "
205 << tpC->size());
206 for (const auto* trk : eleIDtpStored)
207 ATH_MSG_DEBUG("ele ID trk "
208 << trk << " pt = " << trk->pt() * 1e-3
209 << " eta = " << trk->eta() << " phi = " << trk->phi()
210 << " nSi = " << xAOD::EgammaHelpers::numberOfSiHits(trk));
211 for (const auto* trk : eleGSFtpStored)
212 ATH_MSG_DEBUG("ele GSF trk "
213 << trk << " pt = " << trk->pt() * 1e-3
214 << " eta = " << trk->eta() << " phi = " << trk->phi()
215 << " nSi = " << xAOD::EgammaHelpers::numberOfSiHits(trk));
216 for (const xAOD::TrackParticle* trk : *closeByTracks)
217 ATH_MSG_DEBUG("closeby trk "
218 << trk << " pt = " << trk->pt() * 1e-3
219 << " eta = " << trk->eta() << " phi = " << trk->phi()
220 << " nSi = " << xAOD::EgammaHelpers::numberOfSiHits(trk));
221 }
222
223 for (const auto* ele : *eleC) {
224
225 // Electron preselection
226 if (ele->pt() < m_elepTCut ||
227 m_idCut.empty() || !aidCut.isAvailable(*ele) || !aidCut(*ele))
228 continue;
229
230 // Henri's circles
231 if (decorateSimple(dh, closeByTracks, ele, pvtx).isFailure()) {
232 ATH_MSG_ERROR("Cannot decorate the electron with the simple info");
233 return StatusCode::FAILURE;
234 }
235
236 } // loop on electrons to decorate
237
238 return StatusCode::SUCCESS;
239}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
size_type size() const noexcept
Returns the number of elements in the collection.
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_tpContainerName
StatusCode decorateSimple(DecorHandles &dh, std::unique_ptr< ConstDataVector< xAOD::TrackParticleContainer > > &tpC, const xAOD::Electron *ele, const xAOD::Vertex *pvtx) const
SG::ReadHandleKey< xAOD::ElectronContainer > m_containerName
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_tpCName
SG::ReadHandleKey< xAOD::VertexContainer > m_VtxContainerName
const_pointer_type ptr()
Dereference the pointer.
float z0() const
Returns the parameter.
float theta() const
Returns the parameter, which has range 0 to .
virtual FourMom_t p4() const override final
The full 4-momentum of the particle.
@ VIEW_ELEMENTS
this data object is a view, it does not own its elmts
std::size_t numberOfSiHits(const xAOD::TrackParticle *tp)
return the number of Si hits in the track particle
const xAOD::TrackParticle * getOriginalTrackParticle(const xAOD::Electron *el)
Helper function for getting the "Original" Track Particle (i.e before GSF) via the electron.
@ PriVtx
Primary vertex.
ElectronContainer_v1 ElectronContainer
Definition of the current "electron container version".
TrackParticle_v1 TrackParticle
Reference the current persistent version:
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
Electron_v1 Electron
Definition of the current "egamma version".

◆ decorateSimple()

StatusCode DerivationFramework::EGElectronAmbiguityTool::decorateSimple ( DecorHandles & dh,
std::unique_ptr< ConstDataVector< xAOD::TrackParticleContainer > > & tpC,
const xAOD::Electron * ele,
const xAOD::Vertex * pvtx ) const
private

HelixToCircle Main Track Electron

HelixToCircle Other Electron Conv Track

Definition at line 243 of file EGElectronAmbiguityTool.cxx.

248{
249 // This is the GSF electron track
250 const xAOD::TrackParticle* eleGSFtrkP = ele->trackParticle();
251
252 // And the ID one
253 const xAOD::TrackParticle* eleIDtrkP =
255
256 // For the time being, use the ID one, to be consistent when we only find a
257 // good ID to make a conversion and no GSF
258 bool useGSF = false; // hardcoded because it seems we will not use true. Kept
259 // for the time being, as not 100% sure
260 const xAOD::TrackParticle* eletrkP = useGSF ? eleGSFtrkP : eleIDtrkP;
261
262 ATH_MSG_DEBUG("Electron pt = " << ele->pt() * 1e-3 << " eta = " << ele->eta()
263 << " phi = " << ele->phi() << " GSF trk ptr = "
264 << eleGSFtrkP << " ID trk ptr " << eleIDtrkP);
265
266 if (m_isMC) {
267 const xAOD::TruthParticle* truthEl =
269 double tpvr = -1, tpvp = 9e9, tpvz = 9e9;
270 if (truthEl && MC::isElectron(truthEl) &&
271 truthEl->prodVtx() != nullptr) {
272 tpvr = truthEl->prodVtx()->perp();
273 tpvp = truthEl->prodVtx()->phi();
274 tpvz = truthEl->prodVtx()->z();
275 }
276 dh.dtrv(*ele) = tpvr;
277 dh.dtpv(*ele) = tpvp;
278 dh.dtzv(*ele) = tpvz;
279 }
280
281 // Find the closest track particle with opposite charge and a minimum nb of Si
282 // hits
283 const xAOD::TrackParticle* otrkP(nullptr);
284 double detaMin = 9e9;
285 for (const xAOD::TrackParticle* tp : *tpC) {
286 // Keep only opposite charge
287 if (tp->charge() * eletrkP->charge() > 0)
288 continue;
289
290 // Close-by
291 double dR = eletrkP->p4().DeltaR(tp->p4());
292 double dz = std::abs(eletrkP->z0() - tp->z0()) * sin(eletrkP->theta());
293 if (dR >= 0.3 || dz >= m_dzCut)
294 continue;
295
296 double deta = std::abs(eletrkP->eta() - tp->eta());
297 if (deta < detaMin) {
298 otrkP = tp;
299 detaMin = deta;
300 }
301 }
302
303 double rv = -9e9;
304 double pv = -9e9;
305 double mee = -1.;
306 double meeAtVtx = -1.;
307 double sep = -9e9;
308 bool goodConv = false;
309
310 if (otrkP) {
311
312 // To be consistent with the other, use the ID track.
313 TLorentzVector ep4;
314 ep4.SetPtEtaPhiM(eletrkP->pt(), eletrkP->eta(), eletrkP->phi(), ParticleConstants::electronMassInMeV);
315
316 // Maybe could see if a GSF tp exists for this ID tp and use it if yes ?
317 TLorentzVector op4;
318 op4.SetPtEtaPhiM(otrkP->pt(), otrkP->eta(), otrkP->phi(), ParticleConstants::electronMassInMeV);
319
320 // Simple masses
321 mee = (ep4 + op4).M();
322 op4.SetPhi(eletrkP->phi());
323 meeAtVtx = (ep4 + op4).M();
324
325 // And the conversion point
326 std::vector<double> helix1, helix2;
327 helix1.resize(5);
328 helix2.resize(5);
329 helix(eletrkP, pvtx, helix1);
330 helix(otrkP, pvtx, helix2);
331
332 double beta(0.);
333 if (helix1[4] < helix2[4])
334 beta = TMath::PiOver2() - helix1[4];
335 else
336 beta = TMath::PiOver2() - helix2[4];
337
338 double phi1(helix1[4] + beta);
339 if (phi1 > TMath::TwoPi())
340 phi1 -= TMath::TwoPi();
341 if (phi1 < 0.)
342 phi1 += TMath::TwoPi();
343
344 double phi2(helix2[4] + beta);
345 if (phi2 > TMath::TwoPi())
346 phi2 -= TMath::TwoPi();
347 if (phi2 < 0.)
348 phi2 += TMath::TwoPi();
349
351 double r1 = 1 / (2. * std::abs(helix1[1]));
352
353 double charge1(1.);
354 if (helix1[1] < 0.)
355 charge1 = -1.;
356 double rcenter1(helix1[3] / charge1 + r1);
357 double phicenter1(phi1 + TMath::PiOver2() * charge1);
358
359 double x1 = rcenter1 * cos(phicenter1);
360 double y1 = rcenter1 * sin(phicenter1);
361
363 double r2 = 1 / (2. * std::abs(helix2[1]));
364
365 double charge2(1.);
366 if (helix2[1] < 0.)
367 charge2 = -1.;
368 double rcenter2(helix2[3] / charge2 + r2);
369 double phicenter2(phi2 + TMath::PiOver2() * charge2);
370
371 double x2 = rcenter2 * cos(phicenter2);
372 double y2 = rcenter2 * sin(phicenter2);
374
375 double dx(x1 - x2);
376 if (dx < 1e-9 && dx > 0.)
377 dx = 1e-9;
378 if (dx > -1e-9 && dx < 0.)
379 dx = -1e-9;
380 double slope((y1 - y2) / dx);
381 double b(y1 - slope * x1);
382 double alpha(atan(slope));
383 double d(sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2)));
384 // only keeping opposite sign option
385 double separation = d - r1 - r2;
386 double cpx1, cpx2;
387 if (x1 > x2) {
388 cpx1 = x1 - r1 * cos(alpha);
389 cpx2 = x2 + r2 * cos(alpha);
390 } else {
391 cpx1 = x1 + r1 * cos(alpha);
392 cpx2 = x2 - r2 * cos(alpha);
393 }
394
395 double temp1 = (cpx1 + cpx2) / 2;
396 double temp2 = slope * temp1 + b;
397 double convX = cos(beta) * temp1 + sin(beta) * temp2;
398 double convY = -sin(beta) * temp1 + cos(beta) * temp2;
399
400 double dct(helix1[0] - helix2[0]);
401
403 if (std::abs(separation) < m_sepCut && std::abs(dct) < m_dctCut) {
404 goodConv = true;
405 sep = separation;
406 pv = std::atan2(convY, convX);
407 rv = sqrt(convX * convX + convY * convY);
408 if (convX * cos(eletrkP->phi()) + convY * sin(eletrkP->phi()) < 0)
409 rv *= -1.;
410 }
411 } else {
412 dh.dambi(*ele) = -1;
413 }
414 dh.drv(*ele) = rv;
415 dh.dphiv(*ele) = pv;
416 dh.dmee(*ele) = mee;
417 dh.dmeeVtx(*ele) = meeAtVtx;
418 dh.dsep(*ele) = sep;
419 if (goodConv && rv > m_rvECCut && meeAtVtx < m_meeAtVtxECCut)
420 dh.dambi(*ele) = 2;
421 else if (otrkP) {
422 if (mee < m_meeICCut)
423 dh.dambi(*ele) = 1;
424 else
425 dh.dambi(*ele) = 0;
426 }
427 return StatusCode::SUCCESS;
428}
virtual double pt() const override final
The transverse momentum ( ) of the particle.
Definition Egamma_v1.cxx:66
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
Definition Egamma_v1.cxx:71
virtual double phi() const override final
The azimuthal angle ( ) of the particle.
Definition Egamma_v1.cxx:76
const xAOD::TrackParticle * trackParticle(size_t index=0) const
Pointer to the xAOD::TrackParticle/s that match the electron candidate.
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
float charge() const
Returns the charge.
const TruthVertex_v1 * prodVtx() const
The production vertex of this particle.
float z() const
Vertex longitudinal distance along the beam line form the origin.
float phi() const
Vertex azimuthal angle.
float perp() const
Vertex transverse distance from the beam line.
bool isElectron(const T &p)
constexpr double electronMassInMeV
the mass of the electron (in MeV)
const xAOD::TruthParticle * getTruthParticle(const xAOD::IParticle &p)
Return the truthParticle associated to the given IParticle (if any)
TruthParticle_v1 TruthParticle
Typedef to implementation.

◆ initialize()

StatusCode DerivationFramework::EGElectronAmbiguityTool::initialize ( )
finaloverridevirtual

Definition at line 47 of file EGElectronAmbiguityTool.cxx.

48{
49
50 ATH_CHECK(m_containerName.initialize());
51 ATH_CHECK(m_VtxContainerName.initialize());
52 ATH_CHECK(m_tpContainerName.initialize());
53 ATH_CHECK(m_tpCName.initialize());
54
55 ATH_CHECK(m_drv.initialize());
56 ATH_CHECK(m_dphiv.initialize());
57 ATH_CHECK(m_dmee.initialize());
58 ATH_CHECK(m_dmeeVtx.initialize());
59 ATH_CHECK(m_dsep.initialize());
60 ATH_CHECK(m_dambi.initialize());
61 ATH_CHECK(m_dtrv.initialize());
62 ATH_CHECK(m_dtpv.initialize());
63 ATH_CHECK(m_dtzv.initialize());
64
65 return StatusCode::SUCCESS;
66}
#define ATH_CHECK
Evaluate an expression and check for errors.
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_dsep
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_dtzv
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_dambi
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_dmee
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_dtrv
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_dphiv
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_dtpv
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_dmeeVtx
SG::WriteDecorHandleKey< xAOD::ElectronContainer > m_drv

◆ DecorHandles

friend struct DecorHandles
friend

Definition at line 106 of file EGElectronAmbiguityTool.h.

Member Data Documentation

◆ m_containerName

SG::ReadHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_containerName
private
Initial value:
{
this,
"ContainerName",
"Electrons",
"SG key of electron container"
}

Definition at line 40 of file EGElectronAmbiguityTool.h.

40 {
41 this,
42 "ContainerName",
43 "Electrons",
44 "SG key of electron container"
45 };

◆ m_dambi

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_dambi { this, "DFCommonAddAmbiguity", m_containerName, "DFCommonAddAmbiguity", "" }
private

Definition at line 83 of file EGElectronAmbiguityTool.h.

83{ this, "DFCommonAddAmbiguity", m_containerName, "DFCommonAddAmbiguity", "" };

◆ m_dctCut

Gaudi::Property<double> DerivationFramework::EGElectronAmbiguityTool::m_dctCut
private
Initial value:
{
this, "DCTCut", 0.02, "second separation cut"}

Definition at line 125 of file EGElectronAmbiguityTool.h.

125 {
126 this, "DCTCut", 0.02, "second separation cut"};

◆ m_dmee

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_dmee { this, "DFCommonSimpleMee", m_containerName, "DFCommonSimpleMee", "" }
private

Definition at line 74 of file EGElectronAmbiguityTool.h.

74{ this, "DFCommonSimpleMee", m_containerName, "DFCommonSimpleMee", "" };

◆ m_dmeeVtx

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_dmeeVtx { this, "DFCommonSimpleMeeAtVtx", m_containerName, "DFCommonSimpleMeeAtVtx", "" }
private

Definition at line 77 of file EGElectronAmbiguityTool.h.

77{ this, "DFCommonSimpleMeeAtVtx", m_containerName, "DFCommonSimpleMeeAtVtx", "" };

◆ m_dphiv

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_dphiv { this, "DFCommonSimpleConvPhi", m_containerName, "DFCommonSimpleConvPhi", "" }
private

Definition at line 71 of file EGElectronAmbiguityTool.h.

71{ this, "DFCommonSimpleConvPhi", m_containerName, "DFCommonSimpleConvPhi", "" };

◆ m_drv

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_drv { this, "DFCommonSimpleConvRadius", m_containerName, "DFCommonSimpleConvRadius", "" }
private

Definition at line 68 of file EGElectronAmbiguityTool.h.

68{ this, "DFCommonSimpleConvRadius", m_containerName, "DFCommonSimpleConvRadius", "" };

◆ m_dsep

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_dsep { this, "DFCommonSimpleSeparation", m_containerName, "DFCommonSimpleSeparation", "" }
private

Definition at line 80 of file EGElectronAmbiguityTool.h.

80{ this, "DFCommonSimpleSeparation", m_containerName, "DFCommonSimpleSeparation", "" };

◆ m_dtpv

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_dtpv { this, "DFCommonProdTruePhi", m_containerName, "DFCommonProdTruePhi", "" }
private

Definition at line 89 of file EGElectronAmbiguityTool.h.

89{ this, "DFCommonProdTruePhi", m_containerName, "DFCommonProdTruePhi", "" };

◆ m_dtrv

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_dtrv { this, "DFCommonProdTrueRadius", m_containerName, "DFCommonProdTrueRadius", "" }
private

Definition at line 86 of file EGElectronAmbiguityTool.h.

86{ this, "DFCommonProdTrueRadius", m_containerName, "DFCommonProdTrueRadius", "" };

◆ m_dtzv

SG::WriteDecorHandleKey<xAOD::ElectronContainer> DerivationFramework::EGElectronAmbiguityTool::m_dtzv { this, "DFCommonProdTrueZ", m_containerName, "DFCommonProdTrueZ", "" }
private

Definition at line 92 of file EGElectronAmbiguityTool.h.

92{ this, "DFCommonProdTrueZ", m_containerName, "DFCommonProdTrueZ", "" };

◆ m_dzCut

Gaudi::Property<double> DerivationFramework::EGElectronAmbiguityTool::m_dzCut
private
Initial value:
{
this, "dzsinTCut", 0.5, "max dz sinTheta between ele and other tracks"}

Definition at line 129 of file EGElectronAmbiguityTool.h.

129 {
130 this, "dzsinTCut", 0.5, "max dz sinTheta between ele and other tracks"};

◆ m_elepTCut

Gaudi::Property<double> DerivationFramework::EGElectronAmbiguityTool::m_elepTCut
private
Initial value:
{
this, "pTCut", 9000., "minimum pT for an electron to be studied"}

Definition at line 117 of file EGElectronAmbiguityTool.h.

117 {
118 this, "pTCut", 9000., "minimum pT for an electron to be studied"};

◆ m_idCut

Gaudi::Property<std::string> DerivationFramework::EGElectronAmbiguityTool::m_idCut
private
Initial value:
{
this, "idCut", "DFCommonElectronsLHLoose", "minimal quality for an electron to be studied"}

Definition at line 119 of file EGElectronAmbiguityTool.h.

119 {
120 this, "idCut", "DFCommonElectronsLHLoose", "minimal quality for an electron to be studied"};

◆ m_isMC

Gaudi::Property<bool> DerivationFramework::EGElectronAmbiguityTool::m_isMC {this, "isMC", false}
private

Definition at line 114 of file EGElectronAmbiguityTool.h.

114{this, "isMC", false};

◆ m_meeAtVtxECCut

Gaudi::Property<double> DerivationFramework::EGElectronAmbiguityTool::m_meeAtVtxECCut
private
Initial value:
{
this, "meeAtVtxCut", 100, "maximal mass at vertex to be classified as external conversion"}

Definition at line 137 of file EGElectronAmbiguityTool.h.

137 {
138 this, "meeAtVtxCut", 100, "maximal mass at vertex to be classified as external conversion"};

◆ m_meeICCut

Gaudi::Property<double> DerivationFramework::EGElectronAmbiguityTool::m_meeICCut
private
Initial value:
{
this, "meeCut", 100, "maximal mass at primary vertex to be classified as gamma*"}

Definition at line 139 of file EGElectronAmbiguityTool.h.

139 {
140 this, "meeCut", 100, "maximal mass at primary vertex to be classified as gamma*"};

◆ m_nSiCut

Gaudi::Property<unsigned int> DerivationFramework::EGElectronAmbiguityTool::m_nSiCut
private
Initial value:
{
this, "nSiCut", 7, "minimum number of Si hits in the other track"}

Definition at line 123 of file EGElectronAmbiguityTool.h.

123 {
124 this, "nSiCut", 7, "minimum number of Si hits in the other track"};

◆ m_rvECCut

Gaudi::Property<double> DerivationFramework::EGElectronAmbiguityTool::m_rvECCut
private
Initial value:
{
this, "radiusCut", 20, "minimum radius to be classified as external conversion"}

Definition at line 135 of file EGElectronAmbiguityTool.h.

135 {
136 this, "radiusCut", 20, "minimum radius to be classified as external conversion"};

◆ m_sepCut

Gaudi::Property<double> DerivationFramework::EGElectronAmbiguityTool::m_sepCut
private
Initial value:
{
this, "SeparationCut", 1., "first separation cut"}

Definition at line 127 of file EGElectronAmbiguityTool.h.

127 {
128 this, "SeparationCut", 1., "first separation cut"};

◆ m_tpCName

SG::ReadHandleKey<xAOD::TrackParticleContainer> DerivationFramework::EGElectronAmbiguityTool::m_tpCName
private
Initial value:
{
this,
"tpCName",
"GSFTrackParticles",
"SG key of TrackParticleInputContainer"
}

Definition at line 59 of file EGElectronAmbiguityTool.h.

59 {
60 this,
61 "tpCName",
62 "GSFTrackParticles",
63 "SG key of TrackParticleInputContainer"
64 };

◆ m_tpContainerName

SG::ReadHandleKey<xAOD::TrackParticleContainer> DerivationFramework::EGElectronAmbiguityTool::m_tpContainerName
private
Initial value:
{
this,
"tpContainerName",
"InDetTrackParticles",
"SG key of track particles container"
}

Definition at line 52 of file EGElectronAmbiguityTool.h.

52 {
53 this,
54 "tpContainerName",
55 "InDetTrackParticles",
56 "SG key of track particles container"
57 };

◆ m_VtxContainerName

SG::ReadHandleKey<xAOD::VertexContainer> DerivationFramework::EGElectronAmbiguityTool::m_VtxContainerName
private
Initial value:
{
this,
"VtxContainerName",
"PrimaryVertices",
"SG key of vertex container"
}

Definition at line 46 of file EGElectronAmbiguityTool.h.

46 {
47 this,
48 "VtxContainerName",
49 "PrimaryVertices",
50 "SG key of vertex container"
51 };

The documentation for this class was generated from the following files: