CMS 3D CMS Logo

PTDRElectronID.cc
Go to the documentation of this file.
2 
4 
5  // Get all the parameters
6  //baseSetup(conf);
7 
8  quality_ = conf.getParameter<std::string>("electronQuality");
9 
10  useEoverPIn_ = conf.getParameter<std::vector<int> >("useEoverPIn");
11  useDeltaEtaIn_ = conf.getParameter<std::vector<int> >("useDeltaEtaIn");
12  useDeltaPhiIn_ = conf.getParameter<std::vector<int> >("useDeltaPhiIn");
13  useHoverE_ = conf.getParameter<std::vector<int> >("useHoverE");
14  useE9overE25_ = conf.getParameter<std::vector<int> >("useE9overE25");
15  useEoverPOut_ = conf.getParameter<std::vector<int> >("useEoverPOut");
16  useDeltaPhiOut_ = conf.getParameter<std::vector<int> >("useDeltaPhiOut");
17  useInvEMinusInvP_ = conf.getParameter<std::vector<int> >("useInvEMinusInvP");
18  useBremFraction_ = conf.getParameter<std::vector<int> >("useBremFraction");
19  useSigmaEtaEta_ = conf.getParameter<std::vector<int> >("useSigmaEtaEta");
20  useSigmaPhiPhi_ = conf.getParameter<std::vector<int> >("useSigmaPhiPhi");
21  acceptCracks_ = conf.getParameter<std::vector<int> >("acceptCracks");
22 
23  if (quality_=="tight") {
24  cuts_ = conf.getParameter<edm::ParameterSet>("tightEleIDCuts");
25  variables_ = 2 ;
26  } else if (quality_=="medium") {
27  cuts_ = conf.getParameter<edm::ParameterSet>("mediumEleIDCuts");
28  variables_ = 1 ;
29  } else if (quality_=="loose") {
30  cuts_ = conf.getParameter<edm::ParameterSet>("looseEleIDCuts");
31  variables_ = 0 ;
32  } else {
33  throw cms::Exception("Configuration") << "Invalid electronQuality parameter in PTDElectronID: must be tight, medium or loose." ;
34  }
35 }
36 
38  const edm::Event& e ,
39  const edm::EventSetup& es) {
40 
41  //determine which element of the cut arrays in cfi file to read
42  //depending on the electron classification
43  int icut=0;
44  int elClass = electron->classification() ;
45  if (electron->isEB()) //barrel
46  {
47  if (elClass == reco::GsfElectron::GOLDEN) icut=0;
48  if (elClass == reco::GsfElectron::BIGBREM) icut=1;
49  //if (elClass == reco::GsfElectron::NARROW) icut=2;
50  if (elClass == reco::GsfElectron::SHOWERING) icut=3;
51  if (elClass == reco::GsfElectron::GAP) icut=8;
52  }
53  if (electron->isEE()) //endcap
54  {
55  if (elClass == reco::GsfElectron::GOLDEN) icut=4;
56  if (elClass == reco::GsfElectron::BIGBREM) icut=5;
57  //if (elClass == reco::GsfElectron::NARROW) icut=6;
58  if (elClass == reco::GsfElectron::SHOWERING) icut=7;
59  if (elClass == reco::GsfElectron::GAP) icut=8;
60  }
61  if (elClass == reco::GsfElectron::UNKNOWN)
62  {
63  edm::LogError("PTDRElectronID") << "Error: unrecognized electron classification ";
64  return 1.;
65  }
66 
68  if (elClass == reco::GsfElectron::GAP) return 1.;
69 
70  if (useEoverPIn_[variables_]) {
71  double value = electron->eSuperClusterOverP();
72  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPInMax");
73  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPInMin");
74  if (value<mincut[icut] || value>maxcut[icut]) return 0.;
75  }
76 
77  if (useDeltaEtaIn_[variables_]) {
78  double value = electron->deltaEtaSuperClusterTrackAtVtx();
79  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
80  if (fabs(value)>maxcut[icut]) return 0.;
81  }
82 
83  if (useDeltaPhiIn_[variables_]) {
84  double value = electron->deltaPhiSuperClusterTrackAtVtx();
85  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
86  if (fabs(value)>maxcut[icut]) return 0.;
87  }
88 
89  if (useHoverE_[variables_]) {
90  double value = electron->hadronicOverEm();
91  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE");
92  if (value>maxcut[icut]) return 0.;
93  }
94 
95  if (useEoverPOut_[variables_]) {
96  double value = electron->eSeedClusterOverPout();
97  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPOutMax");
98  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPOutMin");
99  if (value<mincut[icut] || value>maxcut[icut]) return 0.;
100  }
101 
102  if (useDeltaPhiOut_[variables_]) {
103  double value = electron->deltaPhiSeedClusterTrackAtCalo();
104  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiOut");
105  if (fabs(value)>maxcut[icut]) return 0.;
106  }
107 
108  if (useInvEMinusInvP_[variables_]) {
109  double value = (1./electron->caloEnergy())-(1./electron->trackMomentumAtVtx().R());
110  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("invEMinusInvP");
111  if (value>maxcut[icut]) return 0.;
112  }
113 
114  if (useBremFraction_[variables_]) {
115  double value = electron->trackMomentumAtVtx().R()-electron->trackMomentumOut().R();
116  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("bremFraction");
117  if (value<mincut[icut]) return 0.;
118  }
119 
120  //EcalClusterLazyTools lazyTools = getClusterShape(e,es);
121  //std::vector<float> vCov = lazyTools.localCovariances(*(electron->superCluster()->seed())) ;
122  //std::vector<float> vCov = lazyTools.covariances(*(electron->superCluster()->seed())) ;
123 
124  if (useE9overE25_[variables_]) {
125  double value = electron->r9()*electron->superCluster()->energy()/electron->e5x5();
126  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("E9overE25");
127  if (fabs(value)<mincut[icut]) return 0.;
128  }
129 
130  if (useSigmaEtaEta_[variables_]) {
131  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMax");
132  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMin");
133  if (electron->sigmaIetaIeta()<mincut[icut] || electron->sigmaIetaIeta()>maxcut[icut]) return 0.;
134  }
135 
136  if (useSigmaPhiPhi_[variables_]) {
137  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMin");
138  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMax");
139  if (electron->sigmaIphiIphi()<mincut[icut] || electron->sigmaIphiIphi()>maxcut[icut]) return 0.;
140  }
141 
142  return 1.;
143 
144 }
float sigmaIphiIphi() const
Definition: GsfElectron.h:441
T getParameter(std::string const &) const
float eSuperClusterOverP() const
Definition: GsfElectron.h:249
edm::ParameterSet cuts_
std::vector< int > useHoverE_
math::XYZVectorF trackMomentumAtVtx() const
Definition: GsfElectron.h:295
std::vector< int > useInvEMinusInvP_
std::vector< int > useDeltaPhiOut_
bool isEE() const
Definition: GsfElectron.h:357
bool isEB() const
Definition: GsfElectron.h:356
std::vector< int > useSigmaPhiPhi_
double result(const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &) override
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:253
float sigmaIetaIeta() const
Definition: GsfElectron.h:440
float hadronicOverEm() const
Definition: GsfElectron.h:495
std::vector< int > useEoverPIn_
float deltaPhiSeedClusterTrackAtCalo() const
Definition: GsfElectron.h:257
math::XYZVectorF trackMomentumOut() const
Definition: GsfElectron.h:297
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:256
float eSeedClusterOverPout() const
Definition: GsfElectron.h:251
Definition: value.py:1
std::vector< int > useDeltaEtaIn_
std::vector< int > useEoverPOut_
std::vector< int > useSigmaEtaEta_
Classification classification() const
Definition: GsfElectron.h:768
float e5x5() const
Definition: GsfElectron.h:444
float r9() const
Definition: GsfElectron.h:445
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:185
std::vector< int > useBremFraction_
std::vector< int > useDeltaPhiIn_
float caloEnergy() const
Definition: GsfElectron.h:863
std::string quality_
std::vector< int > useE9overE25_
void setup(const edm::ParameterSet &conf) override
std::vector< int > acceptCracks_