ATLAS Offline Software
Loading...
Searching...
No Matches
ElectronSelector Class Reference

#include <ElectronSelector.h>

Inheritance diagram for ElectronSelector:
Collaboration diagram for ElectronSelector:

Public Member Functions

 ElectronSelector ()
 ~ElectronSelector ()
 ElectronSelector (const ElectronSelector &)=delete
ElectronSelectoroperator= (const ElectronSelector &)=delete
void setDebug (bool debug)
const xAOD::TrackParticleGetElecNegTrackParticle (size_t i)
const xAOD::TrackParticleGetElecPosTrackParticle (size_t i)
unsigned int GetElectronCollectionSize ()
void Init ()
void PrepareElectronList (const xAOD::ElectronContainer *pxElecContainer)
bool RecordElectron (const xAOD::Electron *)
void SetPtCut (float newpt)

Static Public Member Functions

template<class T>
static CLHEP::Hep3Vector calculateMomentum (const T *pP)
template<class T>
static float EvalInvMass (const T *pxP1, const T *pxP2, float fMass1, float fMass2=invalidAnswer)
template<class T>
static float EvalInvMass (const T *pxP1, const T *pxP2, const T *pxp3, const T *pxP4, float fMass1, float fMass2=-999.9, float fMass3=-999.9, float fMass4=invalidAnswer)
template<class T>
static float EvalDiMuInvMass (const T *pxP1, const T *pxP2)
template<class T>
static float EvalFourMuInvMass (const T *pxP1, const T *pxP2, const T *pxP3, const T *pxP4)
template<class T>
static float EvaluateAngle (const T *pxP1, const T *pxP2)
template<class T>
static float EvalPtDiff (const T *pxP1, const T *pxP2)
template<class T>
static float EvalPhiDiff (const T *pxP1, const T *pxP2)
template<class T>
static float EvalEtaDiff (const T *pxP1, const T *pxP2)
template<class T>
static float EvalPt (const T *pxP1, const T *pxP2)
template<class T>
static float EvalPhi (const T *pxP1, const T *pxP2)
template<class T>
static float EvalEta (const T *pxP1, const T *pxP2)
template<class T>
static float EvalCharge (const T *pxP1, const T *pxP2)
template<class T>
static float EvalTransverseMass (const T *pxP1, float fMETx, float fMETy, float fMass1, float fMass2=invalidAnswer)
template<class T>
static float EvalTransverseMass (const T *pxP1, float fMETx, float fMETy)
template<class T>
static float EvalTransverseMass (const T *pxP1, const T *pxP2, float fMETx, float fMETy, float fMass1, float fMass2=invalidAnswer)
template<class T>
static float EvalTransverseMass (const T *pxP1, const T *pxP2, float fMETx, float fMETy)

Static Public Attributes

static constexpr float invalidAnswer {-999.9f}

Protected Member Functions

virtual void BookHistograms ()

Protected Attributes

unsigned int m_uPassedEvents
std::map< unsigned int, TH1F * > m_x1DHistograms
std::map< unsigned int, TH2F * > m_x2DHistograms
std::map< unsigned int, TProfile * > m_x1DProfHistograms
std::map< unsigned int, TProfile2D * > m_x2DProfHistograms
std::string m_xSampleName

Private Types

typedef EventAnalysis PARENT

Private Member Functions

void Clear ()
bool OrderElectronList ()
bool RetrieveVertices ()
void Register ()

Private Attributes

MsgStream * m_msgStream
const xAOD::Muonm_pxElectron
std::vector< const xAOD::TrackParticle * > m_pxElTrackList
std::vector< const xAOD::TrackParticle * > m_goodElecNegTrackParticleList
std::vector< const xAOD::TrackParticle * > m_goodElecPosTrackParticleList
bool m_doDebug
float m_ptCut
float m_etaCut
AsgElectronLikelihoodToolm_LHTool2015 = nullptr
int m_elecneg1 = 0
int m_elecneg2 = 0
int m_elecpos1 = 0
int m_elecpos2 = 0
const float m_CGeV = 1.0e-3

Static Private Attributes

static std::atomic< unsigned int > s_uNumInstances

Detailed Description

Definition at line 30 of file ElectronSelector.h.

Member Typedef Documentation

◆ PARENT

Definition at line 54 of file ElectronSelector.h.

Constructor & Destructor Documentation

◆ ElectronSelector() [1/2]

ElectronSelector::ElectronSelector ( )

Definition at line 30 of file ElectronSelector.cxx.

30 :
31 m_doDebug ( false ),
32 m_ptCut ( 10. ),
33 m_etaCut ( 2.47 ) // 2.47 is the official acceptance for central electrons. Forward electrons is another story...
34{
36
37 std::stringstream xTmp; xTmp << s_uNumInstances;
38 m_xSampleName = "ElectronSelector_" + xTmp.str();
39
40 m_pxElectron = nullptr;
41
42 m_msgStream = new MsgStream(Athena::getMessageSvc(), "InDetPerformanceMonitoring" );
43}
MsgStream * m_msgStream
const xAOD::Muon * m_pxElectron
static std::atomic< unsigned int > s_uNumInstances
std::string m_xSampleName
IMessageSvc * getMessageSvc(bool quiet=false)

◆ ~ElectronSelector()

ElectronSelector::~ElectronSelector ( )

Definition at line 46 of file ElectronSelector.cxx.

47{
49 delete m_msgStream;
50}

◆ ElectronSelector() [2/2]

ElectronSelector::ElectronSelector ( const ElectronSelector & )
delete

Member Function Documentation

◆ BookHistograms()

void EventAnalysis::BookHistograms ( )
protectedvirtualinherited

Reimplemented in FourMuonEvent, MuonSelector, and ZmumuEvent.

Definition at line 55 of file EventAnalysis.cxx.

56{
57 // This must be overriden by an inheriting class.
58}

◆ calculateMomentum()

template<class T>
CLHEP::Hep3Vector EventAnalysis::calculateMomentum ( const T * pP)
staticinherited

Definition at line 92 of file EventAnalysis.h.

92 {
93 const auto & p4(pP->p4());
94 return CLHEP::Hep3Vector(p4.Px() * EAna::CGeV , p4.Py() * EAna::CGeV, p4.Pz() * EAna::CGeV);
95}
const float CGeV

◆ Clear()

void ElectronSelector::Clear ( )
private

Definition at line 200 of file ElectronSelector.cxx.

201{
202 m_pxElTrackList.clear();
205
206 // -1 means not assigned
207 m_elecneg1 = -1;
208 m_elecneg2 = -1;
209 m_elecpos1 = -1;
210 m_elecpos2 = -1;
211
212 return;
213}
std::vector< const xAOD::TrackParticle * > m_pxElTrackList
std::vector< const xAOD::TrackParticle * > m_goodElecPosTrackParticleList
std::vector< const xAOD::TrackParticle * > m_goodElecNegTrackParticleList

◆ EvalCharge()

template<class T>
float EventAnalysis::EvalCharge ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 239 of file EventAnalysis.h.

240{
241 // Check integrity of inputs.
242 if ( !pxP1 || !pxP2 ) return invalidAnswer;
243 return static_cast<float>( pxP1->charge() + pxP2->charge() );
244}
static constexpr float invalidAnswer

◆ EvalDiMuInvMass()

template<class T>
float EventAnalysis::EvalDiMuInvMass ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 98 of file EventAnalysis.h.

99{
100 // Check integrity of inputs.
101 if ( !pxP1 || !pxP2 ) return invalidAnswer;
102
103 // Evaluate Di-mu invariant mass.
104 return EvalInvMass( pxP1, pxP2, EAna::g_fMuonMass );
105}
static float EvalInvMass(const T *pxP1, const T *pxP2, float fMass1, float fMass2=invalidAnswer)
const float g_fMuonMass

◆ EvalEta()

template<class T>
float EventAnalysis::EvalEta ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 230 of file EventAnalysis.h.

231{
232 // Check integrity of inputs.
233 if ( !pxP1 || !pxP2 ) return invalidAnswer;
234 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
235 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
236 return static_cast<float>( (xTmp1 + xTmp2).pseudoRapidity() );
237}
static CLHEP::Hep3Vector calculateMomentum(const T *pP)

◆ EvalEtaDiff()

template<class T>
float EventAnalysis::EvalEtaDiff ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 202 of file EventAnalysis.h.

203{
204 // Check integrity of inputs.
205 if ( !pxP1 || !pxP2 ) return invalidAnswer;
206 // Evaluate the angle.
207 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
208 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
209 return static_cast<float>( xTmp1.polarAngle(xTmp2) );
210}

◆ EvalFourMuInvMass()

template<class T>
float EventAnalysis::EvalFourMuInvMass ( const T * pxP1,
const T * pxP2,
const T * pxP3,
const T * pxP4 )
staticinherited

Definition at line 124 of file EventAnalysis.h.

125{
126 // Check integrity of inputs.
127 if ( !pxP1 || !pxP2 || !pxP3 || !pxP4) return invalidAnswer;
128
129 // Evaluate invariant mass.
130 return EvalInvMass( pxP1, pxP2, pxP3, pxP4, EAna::g_fMuonMass );
131}

◆ EvalInvMass() [1/2]

template<class T>
float EventAnalysis::EvalInvMass ( const T * pxP1,
const T * pxP2,
const T * pxp3,
const T * pxP4,
float fMass1,
float fMass2 = -999.9,
float fMass3 = -999.9,
float fMass4 = invalidAnswer )
staticinherited

Definition at line 133 of file EventAnalysis.h.

135{
136 // Check integrity of inputs.No tachyons.
137 if ( fMass1 < 0.0f ) return invalidAnswer;
138 if ( !pxP1 || !pxP2 || !pxP3 || !pxP4) return invalidAnswer;
139
140 // Set masses equal if required by user
141 fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
142 fMass3 = ( fMass3 < 0.0f ) ? fMass1 : fMass3;
143 fMass4 = ( fMass4 < 0.0f ) ? fMass1 : fMass4;
144
145 // Evaluate invariant mass.
146 CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, pxP1->p4().Pz() * EAna::CGeV );
147 CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( pxP2->p4().Px() * EAna::CGeV, pxP2->p4().Py() * EAna::CGeV, pxP2->p4().Pz() * EAna::CGeV );
148 CLHEP::Hep3Vector xTmp3 = CLHEP::Hep3Vector( pxP3->p4().Px() * EAna::CGeV, pxP3->p4().Py() * EAna::CGeV, pxP3->p4().Pz() * EAna::CGeV );
149 CLHEP::Hep3Vector xTmp4 = CLHEP::Hep3Vector( pxP4->p4().Px() * EAna::CGeV, pxP4->p4().Py() * EAna::CGeV, pxP4->p4().Pz() * EAna::CGeV );
150
151 CLHEP::HepLorentzVector xLVec; xLVec.setVectM ( xTmp1, fMass1 );
152 CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, fMass2 );
153 CLHEP::HepLorentzVector xLVec3; xLVec3.setVectM( xTmp3, fMass3 );
154 CLHEP::HepLorentzVector xLVec4; xLVec4.setVectM( xTmp4, fMass4 );
155
156 xLVec += xLVec2;
157 xLVec += xLVec3;
158 xLVec += xLVec4;
159
160 return static_cast<float>( xLVec.m() );
161}

◆ EvalInvMass() [2/2]

template<class T>
float EventAnalysis::EvalInvMass ( const T * pxP1,
const T * pxP2,
float fMass1,
float fMass2 = invalidAnswer )
staticinherited

Definition at line 107 of file EventAnalysis.h.

109{
110 // Check integrity of inputs.No tachyons.
111 if ( fMass1 < 0.0f ) return invalidAnswer;
112 if ( !pxP1 || !pxP2 ) return invalidAnswer;
113 // Set masses equal if required by user
114 fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
115 // Evaluate invariant mass.
116 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
117 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
118 CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp1, fMass1 );
119 CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, fMass2 );
120 return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
121}

◆ EvalPhi()

template<class T>
float EventAnalysis::EvalPhi ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 221 of file EventAnalysis.h.

222{
223 // Check integrity of inputs.
224 if ( !pxP1 || !pxP2 ) return invalidAnswer;
225 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
226 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
227 return static_cast<float>( (xTmp1 + xTmp2).phi() );
228}

◆ EvalPhiDiff()

template<class T>
float EventAnalysis::EvalPhiDiff ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 191 of file EventAnalysis.h.

192{
193 // Check integrity of inputs.
194 if ( !pxP1 || !pxP2 ) return invalidAnswer;
195 // Evaluate the angle.
196 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
197 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
198 return static_cast<float>( xTmp1.deltaPhi(xTmp2) );
199}

◆ EvalPt()

template<class T>
float EventAnalysis::EvalPt ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 212 of file EventAnalysis.h.

213{
214 // Check integrity of inputs.
215 if ( !pxP1 || !pxP2 ) return invalidAnswer;
216 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
217 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
218 return static_cast<float>( (xTmp1 + xTmp2).perp() );
219}

◆ EvalPtDiff()

template<class T>
float EventAnalysis::EvalPtDiff ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 174 of file EventAnalysis.h.

175{
176 // Check integrity of inputs.
177 if ( !pxP1 || !pxP2 ) return invalidAnswer;
178 // Evaluate the difference between the momenta. Signed using positive - negative if appropriate.
179 if ( pxP1->charge() > 0.5f )
180 {
181 return static_cast<float>( pxP1->pt() * EAna::CGeV - pxP2->pt() * EAna::CGeV );
182 }
183 else
184 {
185 return static_cast<float>( pxP2->pt() * EAna::CGeV - pxP1->pt() * EAna::CGeV );
186 }
187}

◆ EvalTransverseMass() [1/4]

template<class T>
float EventAnalysis::EvalTransverseMass ( const T * pxP1,
const T * pxP2,
float fMETx,
float fMETy )
staticinherited

Definition at line 272 of file EventAnalysis.h.

273{
274 // Check integrity of inputs.
275 if ( !pxP1 || !pxP2 ) return invalidAnswer;
276 // Evaluate Di-mu invariant mass.
277 return EvalTransverseMass( pxP1, pxP2, fMETx, fMETy, EAna::g_fMuonMass );
278}
static float EvalTransverseMass(const T *pxP1, float fMETx, float fMETy, float fMass1, float fMass2=invalidAnswer)

◆ EvalTransverseMass() [2/4]

template<class T>
float EventAnalysis::EvalTransverseMass ( const T * pxP1,
const T * pxP2,
float fMETx,
float fMETy,
float fMass1,
float fMass2 = invalidAnswer )
staticinherited

Definition at line 280 of file EventAnalysis.h.

282{
283 // Check integrity of inputs.No tachyons.
284 if ( fMass1 < 0.0f ) return invalidAnswer;
285 if ( !pxP1 || !pxP2 ) return invalidAnswer;
286 // Set masses equal if required by user.
287 fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
288 // Evaluate invariant mass.
289 CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, 0.0f );
290 CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( pxP2->p4().Px() * EAna::CGeV, pxP2->p4().Py() * EAna::CGeV, 0.0f );
291 CLHEP::Hep3Vector xTmp12 = xTmp1 + xTmp2;
292 CLHEP::Hep3Vector xTmp3 = CLHEP::Hep3Vector( fMETx, fMETy, 0.0f );
293 CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp12, fMass1 );
294 CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp3, 0.0f );
295 return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
296}

◆ EvalTransverseMass() [3/4]

template<class T>
float EventAnalysis::EvalTransverseMass ( const T * pxP1,
float fMETx,
float fMETy )
staticinherited

Definition at line 247 of file EventAnalysis.h.

248{
249 // Check integrity of inputs.
250 if ( !pxP1 ) return invalidAnswer;
251 // Evaluate Di-mu invariant mass.
252 return EvalInvMass( pxP1, fMETx, fMETy, EAna::g_fMuonMass );
253}

◆ EvalTransverseMass() [4/4]

template<class T>
float EventAnalysis::EvalTransverseMass ( const T * pxP1,
float fMETx,
float fMETy,
float fMass1,
float fMass2 = invalidAnswer )
staticinherited

Definition at line 255 of file EventAnalysis.h.

257{
258 // Check integrity of inputs.No tachyons.
259 if ( fMass1 < 0.0f ) return invalidAnswer;
260 if ( !pxP1 ) return invalidAnswer;
261 // Set masses equal if required by user.
262 fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
263 // Evaluate invariant mass.
264 CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, 0.0f );
265 CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( fMETx, fMETy, 0.0f );
266 CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp1, fMass1 );
267 CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, 0.0f );
268 return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
269}

◆ EvaluateAngle()

template<class T>
float EventAnalysis::EvaluateAngle ( const T * pxP1,
const T * pxP2 )
staticinherited

Definition at line 164 of file EventAnalysis.h.

165{
166 // Check integrity of inputs.
167 if ( !pxP1 || !pxP2 ) return invalidAnswer;
168 // Evaluate the angle.
169 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
170 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
171 return static_cast<float>( xTmp1.angle(xTmp2) );
172}

◆ GetElecNegTrackParticle()

const xAOD::TrackParticle * ElectronSelector::GetElecNegTrackParticle ( size_t i)

Definition at line 349 of file ElectronSelector.cxx.

350{
351 if (i >= m_goodElecNegTrackParticleList.size()) { // requesting out of range electron
352 return nullptr;
353 }
355}

◆ GetElecPosTrackParticle()

const xAOD::TrackParticle * ElectronSelector::GetElecPosTrackParticle ( size_t i)

Definition at line 358 of file ElectronSelector.cxx.

359{
360 if (i >= m_goodElecPosTrackParticleList.size()) { // requesting out of range electron
361 return nullptr;
362 }
364}

◆ GetElectronCollectionSize()

unsigned int ElectronSelector::GetElectronCollectionSize ( )
inline

Definition at line 43 of file ElectronSelector.h.

◆ Init()

void ElectronSelector::Init ( )
virtual

Reimplemented from EventAnalysis.

Definition at line 53 of file ElectronSelector.cxx.

54{
55 (*m_msgStream) << MSG::DEBUG << " -- ElectronSelector::Init -- START -- " << endmsg;
56
57 // PARENT::Init();
58
59 //---Electron Likelihood tool---
60 // m_doIDCuts = true;
61 (*m_msgStream) << MSG::INFO << "ElectronSelector::Init -- Setting up electron LH tool." << endmsg;
62 m_LHTool2015 = new AsgElectronLikelihoodTool ("m_LHTool2015");
63
64 const std::string elecWorkingPoint = "LooseLHElectron"; // "MediumLHElectron" "TightLHElectron"
65
66 if((m_LHTool2015->setProperty("WorkingPoint",elecWorkingPoint.c_str())).isFailure()) {
67 (*m_msgStream) << MSG::WARNING << "Failure loading ConfigFile for electron likelihood tool with working point: " << elecWorkingPoint.c_str() << endmsg;
68 }
69 else {
70 (*m_msgStream) << MSG::INFO << "Loading ConfigFile for electron likelihood tool with working point: " << elecWorkingPoint << ". SUCCESS " << endmsg;
71 }
72
73 // check config files at: https://twiki.cern.ch/twiki/bin/viewauth/AtlasProtected/EGammaIdentificationRun2
74 std::string confDir = "ElectronPhotonSelectorTools/offline/mc20_20210514/ElectronLikelihoodVeryLooseOfflineConfig2017_Smooth.conf";
75 if ( (m_LHTool2015->setProperty("ConfigFile", confDir)).isSuccess()) {
76 (*m_msgStream) << MSG::INFO << "Electron likelihood config ("<< confDir.c_str() << ") setting SUCCESS!" << endmsg;
77 }
78 else {
79 (*m_msgStream) << MSG::WARNING << "Electron likelihood config ("<< confDir.c_str() << ") setting FAILURE" << endmsg;
80 }
81
82 if (m_LHTool2015->initialize().isSuccess()) {
83 (*m_msgStream) << MSG::INFO << "Electron likelihood tool initialize() SUCCESS!" << endmsg;
84 }
85 else {
86 (*m_msgStream) << MSG::WARNING << "Electron likelihood tool initialize() FAILURE!" << endmsg;
87 }
88
89 (*m_msgStream) << MSG::DEBUG << " --ElectronSelector::Init -- COMPLETED -- " << endmsg;
90 return;
91}
#define endmsg
AsgElectronLikelihoodTool * m_LHTool2015

◆ operator=()

ElectronSelector & ElectronSelector::operator= ( const ElectronSelector & )
delete

◆ OrderElectronList()

bool ElectronSelector::OrderElectronList ( )
private

Definition at line 216 of file ElectronSelector.cxx.

217{
218 (*m_msgStream) << MSG::DEBUG << " -- ElectronSelector::OrderElectronList -- START -- list size: " << m_pxElTrackList.size( ) << endmsg;
219
220 bool goodlist = true;
221
222 if (m_pxElTrackList.size() >= 2) { // we need at least 2 electrons
223 double ptMinus1 = 0.;
224 double ptMinus2 = 0.;
225 double ptPlus1 = 0.;
226 double ptPlus2 = 0.;
227
228 int elecnegcount = 0;
229 int elecposcount = 0;
230
231 for (int ielec=0; ielec < (int) m_pxElTrackList.size(); ielec++) {
232 // negative electrons
233 if (m_pxElTrackList.at(ielec)->charge()== -1) { // positive electron
234 elecnegcount++;
235 if (m_pxElTrackList.at(ielec)->pt()> ptMinus1) {
236 // store 1st in 2nd
237 ptMinus2 = ptMinus1;
239 // now store the new one in 1st place
240 ptMinus1 = m_pxElTrackList.at(ielec)->pt();
241 m_elecneg1 = ielec;
242 }
243 else if (m_pxElTrackList.at(ielec)->pt()> ptMinus2) {
244 // store the new one in 2nd place
245 ptMinus2 = m_pxElTrackList.at(ielec)->pt();
246 m_elecneg2 = ielec;
247 }
248 }
249 // positive electrons
250 if (m_pxElTrackList.at(ielec)->charge()== 1) { // positive electron
251 elecposcount++;
252 if (m_pxElTrackList.at(ielec)->pt()> ptPlus1) {
253 // store 1st in 2nd
254 ptPlus2 = ptPlus1;
256 // now store the new one in 1st place
257 ptPlus1 = m_pxElTrackList.at(ielec)->pt();
258 m_elecpos1 = ielec;
259 }
260 else if (m_pxElTrackList.at(ielec)->pt()> ptPlus2) {
261 // store the new one in 2nd place
262 ptPlus2 = m_pxElTrackList.at(ielec)->pt();
263 m_elecpos2 = ielec;
264 }
265 }
266 }
267
268 if (elecposcount == 0 || elecnegcount == 0) {
269 // We need at least one e- and one e+
270 if (m_doDebug) std::cout << " -- ElectronSelector::OrderElectronList -- No opposite charge electrons --> DISCARD ALL ELECTRONS -- " << std::endl;
271 elecposcount = 0;
272 elecnegcount = 0;
273 this->Clear();
274 goodlist = false;
275 }
276
277 if (m_doDebug && elecposcount + elecnegcount >= 2 ){
278 std::cout << " -- ElectronSelector::OrderElectronList -- electron summary list taking " << elecposcount + elecnegcount
279 << " electrons from the input list of " << m_pxElTrackList.size() << " electrons: " << std::endl;
280 if (m_elecneg1 >= 0) std::cout << " leading e-: " << m_elecneg1 << " Pt = " << ptMinus1 << std::endl;
281 if (m_elecneg2 >= 0) std::cout << " second e-: " << m_elecneg2 << " Pt = " << ptMinus2 << std::endl;
282 if (m_elecpos1 >= 0) std::cout << " leading e+: " << m_elecpos1 << " Pt = " << ptPlus1 << std::endl;
283 if (m_elecpos2 >= 0) std::cout << " second e+: " << m_elecpos2 << " Pt = " << ptPlus2 << std::endl;
284 }
285
286 if (elecposcount + elecnegcount >= 2){ // fill the final list of electrons
291 }
292 }
293
294 (*m_msgStream) << MSG::DEBUG << " -- ElectronSelector::OrderElectronList -- COMPLETED -- status: "<< goodlist << std::endl;
295 return goodlist;
296}

◆ PrepareElectronList()

void ElectronSelector::PrepareElectronList ( const xAOD::ElectronContainer * pxElecContainer)

Definition at line 94 of file ElectronSelector.cxx.

95{
96 (*m_msgStream) << MSG::DEBUG << " --ElectronSelector::PrepareElectronList -- START -- " << endmsg;
97 Clear(); // clear current list records
98
99 using electron_iterator = xAOD::ElectronContainer::const_iterator;
100 electron_iterator iter = pxElecContainer->begin();
101 electron_iterator iterEnd = pxElecContainer->end();
102
103 // Loop over the Electrons
104 int electroncount = 0;
105 for(; iter != iterEnd ; ++iter) {
106 electroncount++;
107 (*m_msgStream) << MSG::DEBUG << " -- ElectronSelector::PrepareElectronList -- candiate electron " << electroncount
108 << " has author " << (*iter)->author(xAOD::EgammaParameters::AuthorElectron)
109 << endmsg;
110 const xAOD::Electron * ELE = (*iter);
111 if ( RecordElectron(ELE) ) {
112 (*m_msgStream) << MSG::DEBUG << " -- ElectronSelector::PrepareElectronList -- candiate electron " << electroncount
113 << " is good "
114 << endmsg;
115 }
116 }
117 bool progressingwell = true;
118
119 (*m_msgStream) << MSG::DEBUG << " -- ElectronSelector::PrepareElectronList -- finished recording electrons. "
120 << " recorded electrons: " << m_pxElTrackList.size()
121 << " out of tested electron candidates:" << electroncount << endmsg;
122 if (m_pxElTrackList.size() < 2) progressingwell = false;
123 if (progressingwell) progressingwell = OrderElectronList ();
124 if (progressingwell) progressingwell = RetrieveVertices ();
125
126 if (!progressingwell) {
127 (*m_msgStream) << MSG::DEBUG << " -- ElectronSelector::PrepareElectronList -- FAILED -- this event has not even a good e+e- pair " << endmsg;
128 this->Clear(); // reset the content as it is not going to be used
129 }
130
131 (*m_msgStream) << MSG::DEBUG << " -- ElectronSelector::PrepareElectronList -- COMPLETED -- electroncount -- m_pxElTrackList.size() / all = "
132 << m_pxElTrackList.size() << " / " << electroncount
133 << endmsg;
134 return;
135}
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
bool RecordElectron(const xAOD::Electron *)
const uint16_t AuthorElectron
Object Reconstructed by standard cluster-based algorithm.
Definition EgammaDefs.h:24
Electron_v1 Electron
Definition of the current "egamma version".

◆ RecordElectron()

bool ElectronSelector::RecordElectron ( const xAOD::Electron * thisElec)

Definition at line 138 of file ElectronSelector.cxx.

139{
140 // start assuming electron candidate is good
141 bool electronisgood = true;
142
143 // check the electron satisfies the working point
144 if (!m_LHTool2015->accept(thisElec) ) {
145 electronisgood = false;
146 (*m_msgStream) << MSG::DEBUG << " -- electron fails workingpoint selection -- " << endmsg;
147 }
148
149 //Get the track particle
150 const xAOD::TrackParticle* theTrackParticle = thisElec->trackParticle();
151
152 if (!theTrackParticle) {
153 electronisgood = false;
154 (*m_msgStream) << MSG::DEBUG << " -- electron fails trackparticle -- " << endmsg;
155 }
156
157 if (electronisgood && thisElec->author(xAOD::EgammaParameters::AuthorElectron) != 1) {
158 electronisgood = false;
159 (*m_msgStream) << MSG::DEBUG << " -- electron fails author -- " << thisElec->author(xAOD::EgammaParameters::AuthorElectron) << endmsg;
160 }
161
162 if (electronisgood && theTrackParticle->pt() * m_CGeV < m_ptCut ) { // pt cut given in GeV
163 electronisgood = false;
164 (*m_msgStream) << MSG::DEBUG << " -- electron fails pt cut -- pt= " << theTrackParticle->pt()
165 << " < " << m_ptCut << " (cut value) "
166 << endmsg;
167 }
168
169 const xAOD::CaloCluster* cluster = thisElec->caloCluster();
170 if(!cluster) {
171 electronisgood = false;
172 (*m_msgStream) << MSG::DEBUG << " -- electron candidate has no CaloCluster " << endmsg;
173 }
174
175 if (electronisgood && (cluster->e() * sin(theTrackParticle->theta())) * m_CGeV < m_ptCut) { // cut on et of the cluster
176 electronisgood = false;
177 (*m_msgStream) << MSG::DEBUG << " -- electron fails cluster Et cut -- Et= " << (cluster->e() * cos(theTrackParticle->theta()))* m_CGeV
178 << " < " << m_ptCut << " (cut value) "
179 << endmsg;
180 }
181
182 if (electronisgood && (std::abs(cluster->eta())> m_etaCut || std::abs(theTrackParticle->eta())> m_etaCut) ) { // cut in eta for the cluster and the track
183 electronisgood = false;
184 (*m_msgStream) << MSG::DEBUG << " -- electron fails eta cut -- cluster_eta= " << cluster->eta() << endmsg;
185 }
186
187 if (electronisgood) {
188 // store this electron
189 m_pxElTrackList.push_back(theTrackParticle);
190
191 (*m_msgStream) << MSG::DEBUG << " * RecordElectron * good electron found -> store this electron with pt " << theTrackParticle->pt()
192 << " --> current m_pxElTrackList.size(): " << m_pxElTrackList.size()
193 << std::endl;
194 }
195
196 return electronisgood;
197}
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double e() const
The total energy of the particle.
uint16_t author(uint16_t bitmask=EgammaParameters::AuthorALL) const
Get author.
const xAOD::CaloCluster * caloCluster(size_t index=0) const
Pointer to the xAOD::CaloCluster/s that define the electron candidate.
const xAOD::TrackParticle * trackParticle(size_t index=0) const
Pointer to the xAOD::TrackParticle/s that match the electron candidate.
float theta() const
Returns the parameter, which has range 0 to .
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
TrackParticle_v1 TrackParticle
Reference the current persistent version:

◆ Register()

void EventAnalysis::Register ( )
privateinherited

Definition at line 63 of file EventAnalysis.cxx.

64{
65 ServiceHandle<ITHistSvc> histSvc ("THistSvc", "EventAnalysis");
66
67 // Register histograms in monitoring tool
68 registerHistogramType(*histSvc, m_x1DHistograms, m_xSampleName, "/1dhisto_");
69 registerHistogramType(*histSvc, m_x2DHistograms, m_xSampleName, "/2dhisto_");
70
71 registerHistogramType(*histSvc, m_x1DProfHistograms, m_xSampleName, "/1dprof_");
72 registerHistogramType(*histSvc, m_x2DProfHistograms, m_xSampleName, "/2dprof_");
73}
std::map< unsigned int, TProfile2D * > m_x2DProfHistograms
std::map< unsigned int, TH2F * > m_x2DHistograms
std::map< unsigned int, TH1F * > m_x1DHistograms
std::map< unsigned int, TProfile * > m_x1DProfHistograms

◆ RetrieveVertices()

bool ElectronSelector::RetrieveVertices ( )
private

Definition at line 299 of file ElectronSelector.cxx.

300{
301 if (m_doDebug) std::cout << " -- ElectronSelector::RetrieveVertices -- START -- list size: "
303 << std::endl;
304 bool goodvertices = false;
305 int nverticesfound = 1; // WARNING default must be 0 --> set to 1 for R22 --> needs to be fixed
306
307 if (m_goodElecNegTrackParticleList.size() >= 1 && m_goodElecPosTrackParticleList.size() >= 1) { // we need at least 1 e- and 1 e+
308 // then, check the distances between the e- and e+ vertices, and make sure at least 1 pair comes from same vertex
309 for (size_t ielec = 0; ielec < m_goodElecNegTrackParticleList.size(); ielec++) {
310 // loop on e-
311 // R21 -> R22 SALVA // TrkParticles have no vertex in R22 --> THIS NEEDS A FIX
312 /*
313 if (m_goodElecNegTrackParticleList.at(ielec)->vertex()) {
314 if (m_doDebug) std::cout << " e-(" << ielec <<")->vertex()->v= (" << m_goodElecNegTrackParticleList.at(ielec)->vertex()->x()
315 << ", " << m_goodElecNegTrackParticleList.at(ielec)->vertex()->y()
316 << ", " << m_goodElecNegTrackParticleList.at(ielec)->vertex()->z()
317 << ") " << std::endl;
318
319 // for (unsigned int iposi = 0; iposi < m_goodElecPosTrackParticleList.size(); iposi++) {
320 for (size_t iposi = 0; iposi < m_goodElecPosTrackParticleList.size(); iposi++) {
321 if (m_goodElecPosTrackParticleList.at(iposi)->vertex()) {
322 if (m_doDebug) std::cout << " e+(" << iposi <<")->vertex()->v= (" << m_goodElecPosTrackParticleList.at(iposi)->vertex()->x()
323 << ", " << m_goodElecPosTrackParticleList.at(iposi)->vertex()->y()
324 << ", " << m_goodElecPosTrackParticleList.at(iposi)->vertex()->z()
325 << ") " << std::endl;
326 float delta_x = std::abs( m_goodElecNegTrackParticleList.at(ielec)->vertex()->x()-m_goodElecPosTrackParticleList.at(iposi)->vertex()->x() );
327 float delta_y = std::abs( m_goodElecNegTrackParticleList.at(ielec)->vertex()->y()-m_goodElecPosTrackParticleList.at(iposi)->vertex()->y() );
328 float delta_z = std::abs( m_goodElecNegTrackParticleList.at(ielec)->vertex()->z()-m_goodElecPosTrackParticleList.at(iposi)->vertex()->z() );
329
330 if (delta_x < m_deltaXYcut && delta_y < m_deltaXYcut && delta_z < m_deltaZcut) {
331 nverticesfound++;
332 if (m_doDebug) std::cout << " ELEC-BINGO !!! e+e- pair in same vertex !!! e-[" << ielec
333 << "] e+[" << iposi<< "] count: " << nverticesfound << std::endl;
334 } // vertex is the same
335 } // positron has vertex
336 } // loop on positrons
337 } // electron has vertex
338 */
339 } // loop on electrons (e-)
340 } // at least one e+e- pair
341
342 if (nverticesfound >= 1) goodvertices = true;
343
344 if (m_doDebug) std::cout << " -- ElectronSelector::RetrieveVertices -- COMPLETED -- status: " << goodvertices << std::endl;
345 return goodvertices;
346}

◆ setDebug()

void ElectronSelector::setDebug ( bool debug)
inline

Definition at line 39 of file ElectronSelector.h.

const bool debug

◆ SetPtCut()

void ElectronSelector::SetPtCut ( float newpt)
inline

Definition at line 48 of file ElectronSelector.h.

48{m_ptCut = newpt;}

Member Data Documentation

◆ invalidAnswer

float EventAnalysis::invalidAnswer {-999.9f}
staticconstexprinherited

Definition at line 40 of file EventAnalysis.h.

40{-999.9f};

◆ m_CGeV

const float ElectronSelector::m_CGeV = 1.0e-3
private

Definition at line 88 of file ElectronSelector.h.

◆ m_doDebug

bool ElectronSelector::m_doDebug
private

Definition at line 73 of file ElectronSelector.h.

◆ m_elecneg1

int ElectronSelector::m_elecneg1 = 0
private

Definition at line 82 of file ElectronSelector.h.

◆ m_elecneg2

int ElectronSelector::m_elecneg2 = 0
private

Definition at line 83 of file ElectronSelector.h.

◆ m_elecpos1

int ElectronSelector::m_elecpos1 = 0
private

Definition at line 84 of file ElectronSelector.h.

◆ m_elecpos2

int ElectronSelector::m_elecpos2 = 0
private

Definition at line 85 of file ElectronSelector.h.

◆ m_etaCut

float ElectronSelector::m_etaCut
private

Definition at line 76 of file ElectronSelector.h.

◆ m_goodElecNegTrackParticleList

std::vector<const xAOD::TrackParticle*> ElectronSelector::m_goodElecNegTrackParticleList
private

Definition at line 69 of file ElectronSelector.h.

◆ m_goodElecPosTrackParticleList

std::vector<const xAOD::TrackParticle*> ElectronSelector::m_goodElecPosTrackParticleList
private

Definition at line 70 of file ElectronSelector.h.

◆ m_LHTool2015

AsgElectronLikelihoodTool* ElectronSelector::m_LHTool2015 = nullptr
private

Definition at line 79 of file ElectronSelector.h.

◆ m_msgStream

MsgStream* ElectronSelector::m_msgStream
private

Definition at line 64 of file ElectronSelector.h.

◆ m_ptCut

float ElectronSelector::m_ptCut
private

Definition at line 75 of file ElectronSelector.h.

◆ m_pxElectron

const xAOD::Muon* ElectronSelector::m_pxElectron
private

Definition at line 67 of file ElectronSelector.h.

◆ m_pxElTrackList

std::vector<const xAOD::TrackParticle*> ElectronSelector::m_pxElTrackList
private

Definition at line 68 of file ElectronSelector.h.

◆ m_uPassedEvents

unsigned int EventAnalysis::m_uPassedEvents
protectedinherited

Definition at line 75 of file EventAnalysis.h.

◆ m_x1DHistograms

std::map<unsigned int, TH1F*> EventAnalysis::m_x1DHistograms
protectedinherited

Definition at line 76 of file EventAnalysis.h.

◆ m_x1DProfHistograms

std::map<unsigned int, TProfile*> EventAnalysis::m_x1DProfHistograms
protectedinherited

Definition at line 78 of file EventAnalysis.h.

◆ m_x2DHistograms

std::map<unsigned int, TH2F*> EventAnalysis::m_x2DHistograms
protectedinherited

Definition at line 77 of file EventAnalysis.h.

◆ m_x2DProfHistograms

std::map<unsigned int, TProfile2D*> EventAnalysis::m_x2DProfHistograms
protectedinherited

Definition at line 79 of file EventAnalysis.h.

◆ m_xSampleName

std::string EventAnalysis::m_xSampleName
protectedinherited

Definition at line 81 of file EventAnalysis.h.

◆ s_uNumInstances

std::atomic< unsigned int > ElectronSelector::s_uNumInstances
staticprivate

Definition at line 56 of file ElectronSelector.h.


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