CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronClassification.cc
Go to the documentation of this file.
2 
8 
9 
11 
12 //===================================================================
13 // Author: Federico Ferri - INFN Milano, Bicocca university
14 // 12/2005
15 // See GsfElectron::Classification
16 //===================================================================
17 
18 using namespace reco;
19 
21  classify(electron);
22  // electron.classifyElectron(this);
23  electron.classifyElectron(electronClass_);
24 }
25 
27  {
28  electronClass_ = GsfElectron::UNKNOWN ;
29 
30 
31  reco::SuperClusterRef sclRef=electron.superCluster();
32 
33  // use supercluster energy including f(Ncry) correction
34  float scEnergy=sclRef->energy();
35 
36  // first look whether it's in crack, barrel or endcap
37  if ((!electron.isEB())&&(!electron.isEE()))
38  {
39  edm::LogWarning("") << "ElectronClassification::init(): Undefined electron, eta = " <<
40  electron.eta() << "!!!!" ;
41  return ;
42  }
43 
44  if (electron.isEBEEGap() || electron.isEBEtaGap() || electron.isEERingGap())
45  {
46  electronClass_ = GsfElectron::GAP ;
47  return ;
48  }
49 
50  float pin = electron.trackMomentumAtVtx().R() ;
51  float fbrem = electron.fbrem() ;
52  int nbrem = electron.numberOfBrems() ;
53 
54  // golden
55  if (nbrem == 0 && (pin - scEnergy)/pin < 0.1 && fbrem < 0.5) {
56  electronClass_ = GsfElectron::GOLDEN ;
57  }
58 
59  // big brem
60  else if (nbrem == 0 && (pin - scEnergy)/pin < 0.1 && fbrem > 0.5) {
61  electronClass_ = GsfElectron::BIGBREM ;
62  }
63 
64  // showering
65  else
66  electronClass_ = GsfElectron::SHOWERING ;
67 
68 }
69 
70 /*
71 bool ElectronClassification::isInCrack(float eta) const{
72 
73  return (eta>1.460 && eta<1.558);
74 
75 }
76 
77 bool ElectronClassification::isInEtaGaps(float eta) const{
78 
79  return (eta < 0.018 ||
80  (eta>0.423 && eta<0.461) ||
81  (eta>0.770 && eta<0.806) ||
82  (eta>1.127 && eta<1.163));
83 
84 }
85 
86 bool ElectronClassification::isInPhiGaps(float phi) const{
87 
88  return false;
89 
90 }
91 */
void classifyElectron(Classification myclass)
Definition: GsfElectron.h:533
virtual void correct(reco::GsfElectron &)
bool isEBEtaGap() const
Definition: GsfElectron.h:339
bool isEBEEGap() const
Definition: GsfElectron.h:337
bool isEERingGap() const
Definition: GsfElectron.h:343
float fbrem() const
Definition: GsfElectron.h:528
SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:166
virtual double eta() const
momentum pseudorapidity
bool isEE() const
Definition: GsfElectron.h:335
bool isEB() const
Definition: GsfElectron.h:334
return((rh^lh)&mask)
math::XYZVector trackMomentumAtVtx() const
Definition: GsfElectron.h:244
int numberOfBrems() const
Definition: GsfElectron.h:529
void classify(const reco::GsfElectron &)