CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
ClassBasedElectronID Class Reference

#include <ClassBasedElectronID.h>

Inheritance diagram for ClassBasedElectronID:
ElectronIDAlgo

Public Member Functions

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

Private Attributes

edm::ParameterSet cuts_
 
std::string quality_
 

Additional Inherited Members

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

Detailed Description

Definition at line 6 of file ClassBasedElectronID.h.

Constructor & Destructor Documentation

ClassBasedElectronID::ClassBasedElectronID ( )
inline

Definition at line 10 of file ClassBasedElectronID.h.

10 {};
virtual ClassBasedElectronID::~ClassBasedElectronID ( )
inlinevirtual

Definition at line 12 of file ClassBasedElectronID.h.

12 {};

Member Function Documentation

double ClassBasedElectronID::result ( const reco::GsfElectron electron,
const edm::Event e,
const edm::EventSetup es 
)
virtual

Reimplemented from ElectronIDAlgo.

Definition at line 29 of file ClassBasedElectronID.cc.

References reco::GsfElectron::BIGBREM, reco::LeafCandidate::charge(), reco::GsfElectron::classification(), cuts_, reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), reco::GsfElectron::eSeedClusterOverPout(), reco::GsfElectron::GAP, edm::ParameterSet::getParameter(), reco::GsfElectron::GOLDEN, reco::GsfElectron::hadronicOverEm(), reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), reco::GsfElectron::SHOWERING, reco::GsfElectron::sigmaIetaIeta(), reco::GsfElectron::UNKNOWN, and relativeConstraints::value.

32 {
33 
34  //determine which element of the cut arrays in cfi file to read
35  //depending on the electron classification
36  int icut=0;
37  int elClass = electron->classification() ;
38  if (electron->isEB()) //barrel
39  {
40  if (elClass == reco::GsfElectron::GOLDEN) icut=0;
41  if (elClass == reco::GsfElectron::BIGBREM) icut=1;
42  if (elClass == reco::GsfElectron::SHOWERING) icut=2;
43  if (elClass == reco::GsfElectron::GAP) icut=6;
44  }
45  if (electron->isEE()) //endcap
46  {
47  if (elClass == reco::GsfElectron::GOLDEN) icut=3;
48  if (elClass == reco::GsfElectron::BIGBREM) icut=4;
49  if (elClass == reco::GsfElectron::SHOWERING) icut=5;
50  if (elClass == reco::GsfElectron::GAP) icut=7;
51  }
52  if (elClass == reco::GsfElectron::UNKNOWN)
53  {
54  edm::LogError("ClassBasedElectronID") << "Error: unrecognized electron classification ";
55  return 1.;
56  }
57 
58  bool useDeltaEtaIn = true;
59  bool useSigmaIetaIeta = true;
60  bool useHoverE = true;
61  bool useEoverPOut = true;
62  bool useDeltaPhiInCharge = true;
63 
64  // DeltaEtaIn
65  if (useDeltaEtaIn) {
66  double value = electron->deltaEtaSuperClusterTrackAtVtx();
67  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
68  if (fabs(value)>maxcut[icut]) return 0.;
69  }
70 
71  // SigmaIetaIeta
72  if(useSigmaIetaIeta) {
73  double value = electron->sigmaIetaIeta() ;
74  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaIetaIetaMax");
75  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaIetaIetaMin");
76  if(value<mincut[icut] || value>maxcut[icut]) return 0.;
77  }
78 
79  // H/E
80  if (useHoverE) { //_[variables_]) {
81  double value = electron->hadronicOverEm();
82  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE");
83  if (value>maxcut[icut]) return 0.;
84  } // if use
85 
86  // Eseed/Pout
87  if (useEoverPOut) {
88  double value = electron->eSeedClusterOverPout();
89  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPOutMax");
90  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPOutMin");
91  if (value<mincut[icut] || value>maxcut[icut]) return 0.;
92  }
93 
94  // DeltaPhiIn*Charge
95  if (useDeltaPhiInCharge) {
96  double value1 = electron->deltaPhiSuperClusterTrackAtVtx();
97  double value2 = electron->charge();
98  double value = value1*value2;
99  std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiInChargeMax");
100  std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("deltaPhiInChargeMin");
101  if (value<mincut[icut] || value>maxcut[icut]) return 0.;
102  }
103 
104  return 1.;
105 
106 }
T getParameter(std::string const &) const
bool isEE() const
Definition: GsfElectron.h:351
bool isEB() const
Definition: GsfElectron.h:350
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:247
float sigmaIetaIeta() const
Definition: GsfElectron.h:416
float hadronicOverEm() const
Definition: GsfElectron.h:457
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:250
virtual int charge() const final
electric charge
Definition: LeafCandidate.h:91
float eSeedClusterOverPout() const
Definition: GsfElectron.h:245
Classification classification() const
Definition: GsfElectron.h:694
edm::ParameterSet cuts_
void ClassBasedElectronID::setup ( const edm::ParameterSet conf)
virtual

Reimplemented from ElectronIDAlgo.

Definition at line 4 of file ClassBasedElectronID.cc.

References cuts_, cmsRelvalreport::exit, edm::ParameterSet::getParameter(), quality_, and AlCaHLTBitMon_QueryRunRegistry::string.

6 {
7 
8 
9  // Get all the parameters
10  //baseSetup(conf);
11 
12  quality_ = conf.getParameter<std::string>("electronQuality");
13 
14  if(quality_=="Eff95Cuts") {
15  cuts_ = conf.getParameter<edm::ParameterSet>("Eff95Cuts");
16  }
17 
18  else if(quality_=="Eff90Cuts") {
19  cuts_ = conf.getParameter<edm::ParameterSet>("Eff90Cuts");
20  }
21 
22  else {
23  edm::LogError("ClassBasedElectronID") << "Invalid electronQuality parameter: must be tight, medium or loose." ;
24  exit (1);
25  }
26 
27 } // end of setup
T getParameter(std::string const &) const
edm::ParameterSet cuts_

Member Data Documentation

edm::ParameterSet ClassBasedElectronID::cuts_
private

Definition at line 34 of file ClassBasedElectronID.h.

Referenced by result(), and setup().

std::string ClassBasedElectronID::quality_
private

Definition at line 19 of file ClassBasedElectronID.h.

Referenced by setup().