CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
PTDRElectronID Class Reference

#include <PTDRElectronID.h>

Inheritance diagram for PTDRElectronID:
ElectronIDAlgo

Public Member Functions

 PTDRElectronID ()
 
double result (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &) override
 
void setup (const edm::ParameterSet &conf) override
 
 ~PTDRElectronID () override
 
- Public Member Functions inherited from ElectronIDAlgo
 ElectronIDAlgo ()
 
virtual ~ElectronIDAlgo ()
 

Private Attributes

std::vector< int > acceptCracks_
 
edm::ParameterSet cuts_
 
std::string quality_
 
std::vector< int > useBremFraction_
 
std::vector< int > useDeltaEtaIn_
 
std::vector< int > useDeltaPhiIn_
 
std::vector< int > useDeltaPhiOut_
 
std::vector< int > useE9overE25_
 
std::vector< int > useEoverPIn_
 
std::vector< int > useEoverPOut_
 
std::vector< int > useHoverE_
 
std::vector< int > useInvEMinusInvP_
 
std::vector< int > useSigmaEtaEta_
 
std::vector< int > useSigmaPhiPhi_
 
int variables_
 

Additional Inherited Members

- Protected Attributes inherited from ElectronIDAlgo
edm::InputTag reducedBarrelRecHitCollection_
 
edm::InputTag reducedEndcapRecHitCollection_
 

Detailed Description

Definition at line 6 of file PTDRElectronID.h.

Constructor & Destructor Documentation

◆ PTDRElectronID()

PTDRElectronID::PTDRElectronID ( )
inline

Definition at line 8 of file PTDRElectronID.h.

8 {};

◆ ~PTDRElectronID()

PTDRElectronID::~PTDRElectronID ( )
inlineoverride

Definition at line 10 of file PTDRElectronID.h.

10 {};

Member Function Documentation

◆ result()

double PTDRElectronID::result ( const reco::GsfElectron electron,
const edm::Event e,
const edm::EventSetup es 
)
overridevirtual

Reimplemented from ElectronIDAlgo.

Definition at line 37 of file PTDRElectronID.cc.

37  {
38  //determine which element of the cut arrays in cfi file to read
39  //depending on the electron classification
40  int icut = 0;
41  int elClass = electron->classification();
42  if (electron->isEB()) //barrel
43  {
44  if (elClass == reco::GsfElectron::GOLDEN)
45  icut = 0;
46  if (elClass == reco::GsfElectron::BIGBREM)
47  icut = 1;
48  //if (elClass == reco::GsfElectron::NARROW) icut=2;
49  if (elClass == reco::GsfElectron::SHOWERING)
50  icut = 3;
51  if (elClass == reco::GsfElectron::GAP)
52  icut = 8;
53  }
54  if (electron->isEE()) //endcap
55  {
56  if (elClass == reco::GsfElectron::GOLDEN)
57  icut = 4;
58  if (elClass == reco::GsfElectron::BIGBREM)
59  icut = 5;
60  //if (elClass == reco::GsfElectron::NARROW) icut=6;
61  if (elClass == reco::GsfElectron::SHOWERING)
62  icut = 7;
63  if (elClass == reco::GsfElectron::GAP)
64  icut = 8;
65  }
66  if (elClass == reco::GsfElectron::UNKNOWN) {
67  edm::LogError("PTDRElectronID") << "Error: unrecognized electron classification ";
68  return 1.;
69  }
70 
72  if (elClass == reco::GsfElectron::GAP)
73  return 1.;
74 
75  if (useEoverPIn_[variables_]) {
76  double value = electron->eSuperClusterOverP();
77  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPInMax");
78  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPInMin");
79  if (value < mincut[icut] || value > maxcut[icut])
80  return 0.;
81  }
82 
84  double value = electron->deltaEtaSuperClusterTrackAtVtx();
85  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
86  if (fabs(value) > maxcut[icut])
87  return 0.;
88  }
89 
91  double value = electron->deltaPhiSuperClusterTrackAtVtx();
92  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
93  if (fabs(value) > maxcut[icut])
94  return 0.;
95  }
96 
97  if (useHoverE_[variables_]) {
98  double value = electron->hadronicOverEm();
99  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE");
100  if (value > maxcut[icut])
101  return 0.;
102  }
103 
104  if (useEoverPOut_[variables_]) {
105  double value = electron->eSeedClusterOverPout();
106  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPOutMax");
107  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPOutMin");
108  if (value < mincut[icut] || value > maxcut[icut])
109  return 0.;
110  }
111 
113  double value = electron->deltaPhiSeedClusterTrackAtCalo();
114  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiOut");
115  if (fabs(value) > maxcut[icut])
116  return 0.;
117  }
118 
120  double value = (1. / electron->caloEnergy()) - (1. / electron->trackMomentumAtVtx().R());
121  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("invEMinusInvP");
122  if (value > maxcut[icut])
123  return 0.;
124  }
125 
127  double value = electron->trackMomentumAtVtx().R() - electron->trackMomentumOut().R();
128  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("bremFraction");
129  if (value < mincut[icut])
130  return 0.;
131  }
132 
133  //EcalClusterLazyTools lazyTools = getClusterShape(e,es);
134  //std::vector<float> vCov = lazyTools.localCovariances(*(electron->superCluster()->seed())) ;
135  //std::vector<float> vCov = lazyTools.covariances(*(electron->superCluster()->seed())) ;
136 
137  if (useE9overE25_[variables_]) {
138  double value = electron->r9() * electron->superCluster()->energy() / electron->e5x5();
139  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("E9overE25");
140  if (fabs(value) < mincut[icut])
141  return 0.;
142  }
143 
145  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMax");
146  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMin");
147  if (electron->sigmaIetaIeta() < mincut[icut] || electron->sigmaIetaIeta() > maxcut[icut])
148  return 0.;
149  }
150 
152  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMin");
153  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMax");
154  if (electron->sigmaIphiIphi() < mincut[icut] || electron->sigmaIphiIphi() > maxcut[icut])
155  return 0.;
156  }
157 
158  return 1.;
159 }

References acceptCracks_, reco::GsfElectron::BIGBREM, cuts_, metsig::electron, reco::GsfElectron::GAP, edm::ParameterSet::getParameter(), reco::GsfElectron::GOLDEN, reco::GsfElectron::SHOWERING, reco::GsfElectron::UNKNOWN, useBremFraction_, useDeltaEtaIn_, useDeltaPhiIn_, useDeltaPhiOut_, useE9overE25_, useEoverPIn_, useEoverPOut_, useHoverE_, useInvEMinusInvP_, useSigmaEtaEta_, useSigmaPhiPhi_, and variables_.

◆ setup()

void PTDRElectronID::setup ( const edm::ParameterSet conf)
overridevirtual

Reimplemented from ElectronIDAlgo.

Definition at line 3 of file PTDRElectronID.cc.

3  {
4  // Get all the parameters
5  //baseSetup(conf);
6 
7  quality_ = conf.getParameter<std::string>("electronQuality");
8 
9  useEoverPIn_ = conf.getParameter<std::vector<int> >("useEoverPIn");
10  useDeltaEtaIn_ = conf.getParameter<std::vector<int> >("useDeltaEtaIn");
11  useDeltaPhiIn_ = conf.getParameter<std::vector<int> >("useDeltaPhiIn");
12  useHoverE_ = conf.getParameter<std::vector<int> >("useHoverE");
13  useE9overE25_ = conf.getParameter<std::vector<int> >("useE9overE25");
14  useEoverPOut_ = conf.getParameter<std::vector<int> >("useEoverPOut");
15  useDeltaPhiOut_ = conf.getParameter<std::vector<int> >("useDeltaPhiOut");
16  useInvEMinusInvP_ = conf.getParameter<std::vector<int> >("useInvEMinusInvP");
17  useBremFraction_ = conf.getParameter<std::vector<int> >("useBremFraction");
18  useSigmaEtaEta_ = conf.getParameter<std::vector<int> >("useSigmaEtaEta");
19  useSigmaPhiPhi_ = conf.getParameter<std::vector<int> >("useSigmaPhiPhi");
20  acceptCracks_ = conf.getParameter<std::vector<int> >("acceptCracks");
21 
22  if (quality_ == "tight") {
23  cuts_ = conf.getParameter<edm::ParameterSet>("tightEleIDCuts");
24  variables_ = 2;
25  } else if (quality_ == "medium") {
26  cuts_ = conf.getParameter<edm::ParameterSet>("mediumEleIDCuts");
27  variables_ = 1;
28  } else if (quality_ == "loose") {
29  cuts_ = conf.getParameter<edm::ParameterSet>("looseEleIDCuts");
30  variables_ = 0;
31  } else {
32  throw cms::Exception("Configuration")
33  << "Invalid electronQuality parameter in PTDElectronID: must be tight, medium or loose.";
34  }
35 }

References acceptCracks_, cuts_, Exception, edm::ParameterSet::getParameter(), quality_, AlCaHLTBitMon_QueryRunRegistry::string, useBremFraction_, useDeltaEtaIn_, useDeltaPhiIn_, useDeltaPhiOut_, useE9overE25_, useEoverPIn_, useEoverPOut_, useHoverE_, useInvEMinusInvP_, useSigmaEtaEta_, useSigmaPhiPhi_, and variables_.

Member Data Documentation

◆ acceptCracks_

std::vector<int> PTDRElectronID::acceptCracks_
private

Definition at line 29 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ cuts_

edm::ParameterSet PTDRElectronID::cuts_
private

Definition at line 31 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ quality_

std::string PTDRElectronID::quality_
private

Definition at line 16 of file PTDRElectronID.h.

Referenced by setup().

◆ useBremFraction_

std::vector<int> PTDRElectronID::useBremFraction_
private

Definition at line 26 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useDeltaEtaIn_

std::vector<int> PTDRElectronID::useDeltaEtaIn_
private

Definition at line 19 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useDeltaPhiIn_

std::vector<int> PTDRElectronID::useDeltaPhiIn_
private

Definition at line 20 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useDeltaPhiOut_

std::vector<int> PTDRElectronID::useDeltaPhiOut_
private

Definition at line 24 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useE9overE25_

std::vector<int> PTDRElectronID::useE9overE25_
private

Definition at line 22 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useEoverPIn_

std::vector<int> PTDRElectronID::useEoverPIn_
private

Definition at line 18 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useEoverPOut_

std::vector<int> PTDRElectronID::useEoverPOut_
private

Definition at line 23 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useHoverE_

std::vector<int> PTDRElectronID::useHoverE_
private

Definition at line 21 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useInvEMinusInvP_

std::vector<int> PTDRElectronID::useInvEMinusInvP_
private

Definition at line 25 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useSigmaEtaEta_

std::vector<int> PTDRElectronID::useSigmaEtaEta_
private

Definition at line 27 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ useSigmaPhiPhi_

std::vector<int> PTDRElectronID::useSigmaPhiPhi_
private

Definition at line 28 of file PTDRElectronID.h.

Referenced by result(), and setup().

◆ variables_

int PTDRElectronID::variables_
private

Definition at line 33 of file PTDRElectronID.h.

Referenced by result(), and setup().

PTDRElectronID::variables_
int variables_
Definition: PTDRElectronID.h:33
PTDRElectronID::useE9overE25_
std::vector< int > useE9overE25_
Definition: PTDRElectronID.h:22
PTDRElectronID::useHoverE_
std::vector< int > useHoverE_
Definition: PTDRElectronID.h:21
PTDRElectronID::acceptCracks_
std::vector< int > acceptCracks_
Definition: PTDRElectronID.h:29
PTDRElectronID::useSigmaPhiPhi_
std::vector< int > useSigmaPhiPhi_
Definition: PTDRElectronID.h:28
PTDRElectronID::useDeltaPhiOut_
std::vector< int > useDeltaPhiOut_
Definition: PTDRElectronID.h:24
PTDRElectronID::useEoverPIn_
std::vector< int > useEoverPIn_
Definition: PTDRElectronID.h:18
reco::GsfElectron::SHOWERING
Definition: GsfElectron.h:724
metsig::electron
Definition: SignAlgoResolutions.h:48
reco::GsfElectron::GOLDEN
Definition: GsfElectron.h:724
PTDRElectronID::useDeltaEtaIn_
std::vector< int > useDeltaEtaIn_
Definition: PTDRElectronID.h:19
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
reco::GsfElectron::UNKNOWN
Definition: GsfElectron.h:724
edm::ParameterSet
Definition: ParameterSet.h:36
edm::LogError
Definition: MessageLogger.h:183
PTDRElectronID::useSigmaEtaEta_
std::vector< int > useSigmaEtaEta_
Definition: PTDRElectronID.h:27
PTDRElectronID::useEoverPOut_
std::vector< int > useEoverPOut_
Definition: PTDRElectronID.h:23
value
Definition: value.py:1
PTDRElectronID::useBremFraction_
std::vector< int > useBremFraction_
Definition: PTDRElectronID.h:26
reco::GsfElectron::BIGBREM
Definition: GsfElectron.h:724
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Exception
Definition: hltDiff.cc:246
reco::GsfElectron::GAP
Definition: GsfElectron.h:724
PTDRElectronID::useDeltaPhiIn_
std::vector< int > useDeltaPhiIn_
Definition: PTDRElectronID.h:20
PTDRElectronID::cuts_
edm::ParameterSet cuts_
Definition: PTDRElectronID.h:31
PTDRElectronID::quality_
std::string quality_
Definition: PTDRElectronID.h:16
PTDRElectronID::useInvEMinusInvP_
std::vector< int > useInvEMinusInvP_
Definition: PTDRElectronID.h:25