Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
src
RecoEgamma
EgammaElectronAlgos
src
ElectronClassification.cc
Go to the documentation of this file.
1
#include "
RecoEgamma/EgammaElectronAlgos/interface/ElectronClassification.h
"
2
3
#include "
DataFormats/EgammaReco/interface/SuperCluster.h
"
4
5
6
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
7
8
//===================================================================
9
// Author: Federico Ferri - INFN Milano, Bicocca university
10
// 12/2005
11
// See GsfElectron::Classification
12
//===================================================================
13
14
using namespace
reco
;
15
16
void
ElectronClassification::classify
(
GsfElectron
&
electron
)
17
{
18
if
((!electron.
isEB
())&&(!electron.
isEE
()))
19
{
20
edm::LogWarning
(
""
)
21
<<
"ElectronClassification::init(): Undefined electron, eta = "
22
<< electron.
eta
() <<
"!!!!"
;
23
electron.
setClassification
(
GsfElectron::UNKNOWN
) ;
24
return
;
25
}
26
27
if
( electron.
isEBEEGap
() || electron.
isEBEtaGap
() || electron.
isEERingGap
() )
28
{
29
electron.
setClassification
(
GsfElectron::GAP
) ;
30
return
;
31
}
32
33
//float pin = electron.trackMomentumAtVtx().R() ;
34
float
fbrem = electron.
trackFbrem
() ;
35
int
nbrem = electron.
numberOfBrems
() ;
36
37
if
(nbrem == 0 && fbrem < 0.5)
// part (pin - scEnergy)/pin < 0.1 removed - M.D.
38
{ electron.
setClassification
(
GsfElectron::GOLDEN
) ; }
39
else
if
(nbrem == 0 && fbrem >= 0.5)
// part (pin - scEnergy)/pin < 0.1 removed - M.D.
40
{ electron.
setClassification
(
GsfElectron::BIGBREM
) ; }
41
else
42
{ electron.
setClassification
(
GsfElectron::SHOWERING
) ; }
43
44
}
45
46
void
ElectronClassification::refineWithPflow
(
GsfElectron
&
electron
)
47
{
48
if
((!electron.
isEB
())&&(!electron.
isEE
()))
49
{
return
; }
50
51
if
( electron.
isEBEEGap
() || electron.
isEBEtaGap
() || electron.
isEERingGap
() )
52
{
return
; }
53
54
if
((electron.
pfSuperClusterFbrem
()-electron.
trackFbrem
())>=0.15)
55
{ electron.
setClassification
(
GsfElectron::BADTRACK
) ; }
56
}
57
reco::GsfElectron::UNKNOWN
Definition:
GsfElectron.h:673
MessageLogger.h
reco::GsfElectron::GOLDEN
Definition:
GsfElectron.h:673
reco::GsfElectron::BIGBREM
Definition:
GsfElectron.h:673
reco::GsfElectron::isEBEtaGap
bool isEBEtaGap() const
Definition:
GsfElectron.h:352
reco::GsfElectron::isEBEEGap
bool isEBEEGap() const
Definition:
GsfElectron.h:350
reco::GsfElectron
Definition:
GsfElectron.h:37
edm::LogWarning
Definition:
MessageLogger.h:140
reco::GsfElectron::isEERingGap
bool isEERingGap() const
Definition:
GsfElectron.h:356
reco::GsfElectron::BADTRACK
Definition:
GsfElectron.h:673
dt_dqm_sourceclient_common_cff.reco
tuple reco
Definition:
dt_dqm_sourceclient_common_cff.py:105
reco::GsfElectron::GAP
Definition:
GsfElectron.h:673
reco::GsfElectron::isEE
bool isEE() const
Definition:
GsfElectron.h:348
reco::GsfElectron::isEB
bool isEB() const
Definition:
GsfElectron.h:347
reco::GsfElectron::trackFbrem
float trackFbrem() const
Definition:
GsfElectron.h:676
reco::GsfElectron::numberOfBrems
int numberOfBrems() const
Definition:
GsfElectron.h:683
reco::LeafCandidate::eta
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
Definition:
LeafCandidate.h:184
ElectronClassification::classify
void classify(reco::GsfElectron &)
Definition:
ElectronClassification.cc:16
metsig::electron
Definition:
SignAlgoResolutions.h:40
reco::GsfElectron::setClassification
void setClassification(Classification myclass)
Definition:
GsfElectron.h:691
ElectronClassification.h
reco::return
return(e1-e2)*(e1-e2)+dp *dp
ElectronClassification::refineWithPflow
void refineWithPflow(reco::GsfElectron &)
Definition:
ElectronClassification.cc:46
reco::GsfElectron::pfSuperClusterFbrem
float pfSuperClusterFbrem() const
Definition:
GsfElectron.h:678
reco::GsfElectron::SHOWERING
Definition:
GsfElectron.h:673
SuperCluster.h
Generated for CMSSW Reference Manual by
1.8.5