CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFElectronSelector.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PatUtils_interface_PFElectronSelector_h
2 #define PhysicsTools_PatUtils_interface_PFElectronSelector_h
3 
4 #ifndef __GCCXML__
6 #endif
10 
13 
14 #include <iostream>
15 
16 class PFElectronSelector : public Selector<pat::Electron> {
17 
18  public: // interface
19 
20  bool verbose_;
21 
23 
25 
26 #ifndef __GCCXML__
28  PFElectronSelector( parameters )
29  {}
30 #endif
31 
33 
34  verbose_ = false;
35 
36  std::string versionStr = parameters.getParameter<std::string>("version");
37 
39 
40  if ( versionStr == "SPRING11" ) {
41  version = SPRING11;
42  }
43  else {
44  throw cms::Exception("InvalidInput") << "Expect version to be one of SPRING11" << std::endl;
45  }
46 
47 
48  initialize( version,
49  parameters.getParameter<double>("MVA"),
50  parameters.getParameter<double>("D0") ,
51  parameters.getParameter<int> ("MaxMissingHits"),
52  parameters.getParameter<std::string> ("electronIDused"),
53  parameters.getParameter<bool> ("ConversionRejection"),
54  parameters.getParameter<double>("PFIso")
55  );
56  if ( parameters.exists("cutsToIgnore") )
57  setIgnoredCuts( parameters.getParameter<std::vector<std::string> >("cutsToIgnore") );
58 
60 
61  }
62 
64  double mva = 0.4,
65  double d0 = 0.02,
66  int nMissingHits = 1,
67  std::string eidUsed = "eidTightMC",
68  bool convRej = true,
69  double pfiso = 0.15 )
70  {
71  version_ = version;
72 
73  // size_t found;
74  // found = eidUsed.find("NONE");
75  // if ( found != string::npos)
76  electronIDvalue_ = eidUsed;
77 
78  push_back("D0", d0 );
79  push_back("MaxMissingHits", nMissingHits );
80  push_back("electronID");
81  push_back("ConversionRejection" );
82  push_back("PFIso", pfiso );
83  push_back("MVA", mva );
84 
85  set("D0");
86  set("MaxMissingHits");
87  set("electronID");
88  set("ConversionRejection", convRej);
89  set("PFIso");
90  set("MVA");
91 
92  indexD0_ = index_type(&bits_, "D0" );
93  indexMaxMissingHits_ = index_type(&bits_, "MaxMissingHits" );
94  indexElectronId_ = index_type(&bits_, "electronID" );
95  indexConvRej_ = index_type(&bits_, "ConversionRejection" );
96  indexPFIso_ = index_type(&bits_, "PFIso" );
97  indexMVA_ = index_type(&bits_, "MVA" );
98  }
99 
100  // Allow for multiple definitions of the cuts.
102  {
103  if (version_ == SPRING11 ) return spring11Cuts(electron, ret);
104  else {
105  return false;
106  }
107  }
108 
110 
111  // cuts based on top group L+J synchronization exercise
113  {
114 
115  ret.set(false);
116 
117  double mva = electron.mva_e_pi();
118  double missingHits = electron.gsfTrack()->hitPattern().numberOfHits(reco::HitPattern::MISSING_INNER_HITS);
119  double corr_d0 = electron.dB();
120 
121  // in >= 39x conversion rejection variables are accessible from Gsf electron
122  Double_t dist = electron.convDist(); // default value is -9999 if conversion partner not found
123  Double_t dcot = electron.convDcot(); // default value is -9999 if conversion partner not found
124  bool isNotConv = !(fabs(dist) < 0.02 && fabs(dcot) < 0.02);
125 
126  int bitWiseResults = (int) electron.electronID( electronIDvalue_ );
127  bool electronIDboolean = ((bitWiseResults & 1) == 1 );
128 
129  double chIso = electron.userIsolation(pat::PfChargedHadronIso);
130  double nhIso = electron.userIsolation(pat::PfNeutralHadronIso);
131  double gIso = electron.userIsolation(pat::PfGammaIso);
132  double et = electron.et() ;
133 
134  double pfIso = (chIso + nhIso + gIso) / et;
135 
136  if ( missingHits <= cut(indexMaxMissingHits_, double()) || ignoreCut(indexMaxMissingHits_) ) passCut(ret, indexMaxMissingHits_ );
137  if ( fabs(corr_d0) < cut(indexD0_, double()) || ignoreCut(indexD0_) ) passCut(ret, indexD0_ );
138  if ( isNotConv || ignoreCut(indexConvRej_) ) passCut(ret, indexConvRej_ );
139  if ( pfIso < cut(indexPFIso_, double()) || ignoreCut(indexPFIso_) ) passCut(ret, indexPFIso_ );
140  if ( mva > cut(indexMVA_, double()) || ignoreCut(indexMVA_) ) passCut(ret, indexMVA_ );
141  if ( electronIDboolean || ignoreCut(indexElectronId_) ) passCut(ret, indexElectronId_);
142  setIgnored(ret);
143  return (bool)ret;
144  }
145 
146  private: // member variables
147 
149 
157 
159 };
160 
161 #endif
index_type indexMaxMissingHits_
T getParameter(std::string const &) const
void set(std::string const &s, bool val=true)
Set a given selection cut, on or off.
Definition: Selector.h:105
dictionary parameters
Definition: Parameters.py:2
bool spring11Cuts(const pat::Electron &electron, pat::strbitset &ret)
virtual double et() const
transverse energy
std::string electronIDvalue_
bool exists(std::string const &parameterName) const
checks if a parameter exists
pat::strbitset::index_type index_type
Definition: Selector.h:29
void setIgnored(pat::strbitset &ret)
set ignored bits
Definition: Selector.h:224
pat::strbitset retInternal_
internal ret if users don&#39;t care about return bits
Definition: Selector.h:287
float convDist() const
Definition: GsfElectron.h:553
float convDcot() const
Definition: GsfElectron.h:554
pat::strbitset bits_
the bitset indexed by strings
Definition: Selector.h:286
void passCut(pat::strbitset &ret, std::string const &s)
Passing cuts.
Definition: Selector.h:176
reco::GsfTrackRef gsfTrack() const
override the reco::GsfElectron::gsfTrack method, to access the internal storage of the supercluster ...
bool ignoreCut(std::string const &s) const
ignore the cut at index &quot;s&quot;
Definition: Selector.h:159
float electronID(const std::string &name) const
Returns a specific electron ID associated to the pat::Electron given its name.
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:46
PFElectronSelector(edm::ParameterSet const &parameters)
Functor that operates on &lt;T&gt;
Definition: Selector.h:24
float userIsolation(IsolationKeys key) const
Definition: Lepton.h:49
PFElectronSelector(edm::ParameterSet const &parameters, edm::ConsumesCollector &&iC)
bool operator()(const pat::Electron &electron, pat::strbitset &ret)
This provides the interface for base classes to select objects.
strbitset & set(bool val=true)
set method of all bits
Definition: strbitset.h:144
Analysis-level electron class.
Definition: Electron.h:52
float mva_e_pi() const
Definition: GsfElectron.h:629
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:212
double dB(IPTYPE type) const
Impact parameter wrt primary vertex or beamspot.
void setIgnoredCuts(std::vector< std::string > const &bitsToIgnore)
set the bits to ignore from a vector
Definition: Selector.h:167
void initialize(Version_t version, double mva=0.4, double d0=0.02, int nMissingHits=1, std::string eidUsed="eidTightMC", bool convRej=true, double pfiso=0.15)
int cut(index_type const &i, int val) const
Access the int cut values at index &quot;s&quot;.
Definition: Selector.h:194