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 
7 
10 
11 #include <iostream>
12 
13 class PFElectronSelector : public Selector<pat::Electron> {
14 
15  public: // interface
16 
17  bool verbose_;
18 
20 
22 
24 
25  verbose_ = false;
26 
27  std::string versionStr = parameters.getParameter<std::string>("version");
28 
30 
31  if ( versionStr == "SPRING11" ) {
32  version = SPRING11;
33  }
34  else {
35  throw cms::Exception("InvalidInput") << "Expect version to be one of SPRING11" << std::endl;
36  }
37 
38 
39  initialize( version,
40  parameters.getParameter<double>("MVA"),
41  parameters.getParameter<double>("D0") ,
42  parameters.getParameter<int> ("MaxMissingHits"),
43  parameters.getParameter<std::string> ("electronIDused"),
44  parameters.getParameter<bool> ("ConversionRejection"),
45  parameters.getParameter<double>("PFIso")
46  );
47  if ( parameters.exists("cutsToIgnore") )
48  setIgnoredCuts( parameters.getParameter<std::vector<std::string> >("cutsToIgnore") );
49 
51 
52  }
53 
55  double mva = 0.4,
56  double d0 = 0.02,
57  int nMissingHits = 1,
58  std::string eidUsed = "eidTightMC",
59  bool convRej = true,
60  double pfiso = 0.15 )
61  {
62  version_ = version;
63 
64  // size_t found;
65  // found = eidUsed.find("NONE");
66  // if ( found != string::npos)
67  electronIDvalue_ = eidUsed;
68 
69  push_back("D0", d0 );
70  push_back("MaxMissingHits", nMissingHits );
71  push_back("electronID");
72  push_back("ConversionRejection" );
73  push_back("PFIso", pfiso );
74  push_back("MVA", mva );
75 
76  set("D0");
77  set("MaxMissingHits");
78  set("electronID");
79  set("ConversionRejection", convRej);
80  set("PFIso");
81  set("MVA");
82 
83  indexD0_ = index_type(&bits_, "D0" );
84  indexMaxMissingHits_ = index_type(&bits_, "MaxMissingHits" );
85  indexElectronId_ = index_type(&bits_, "electronID" );
86  indexConvRej_ = index_type(&bits_, "ConversionRejection" );
87  indexPFIso_ = index_type(&bits_, "PFIso" );
88  indexMVA_ = index_type(&bits_, "MVA" );
89  }
90 
91  // Allow for multiple definitions of the cuts.
93  {
94  if (version_ == SPRING11 ) return spring11Cuts(electron, ret);
95  else {
96  return false;
97  }
98  }
99 
101 
102  // cuts based on top group L+J synchronization exercise
104  {
105 
106  ret.set(false);
107 
108  double mva = electron.mva();
109  double missingHits = electron.gsfTrack()->trackerExpectedHitsInner().numberOfHits() ;
110  double corr_d0 = electron.dB();
111 
112  // in >= 39x conversion rejection variables are accessible from Gsf electron
113  Double_t dist = electron.convDist(); // default value is -9999 if conversion partner not found
114  Double_t dcot = electron.convDcot(); // default value is -9999 if conversion partner not found
115  bool isConv = fabs(dist) < 0.02 && fabs(dcot) < 0.02;
116 
117  int bitWiseResults = (int) electron.electronID( electronIDvalue_ );
118  bool electronIDboolean = ((bitWiseResults & 1) == 1 );
119 
120  double chIso = electron.userIsolation(pat::PfChargedHadronIso);
121  double nhIso = electron.userIsolation(pat::PfNeutralHadronIso);
122  double gIso = electron.userIsolation(pat::PfGammaIso);
123  double et = electron.et() ;
124 
125  double pfIso = (chIso + nhIso + gIso) / et;
126 
127  if ( missingHits <= cut(indexMaxMissingHits_, double()) || ignoreCut(indexMaxMissingHits_) ) passCut(ret, indexMaxMissingHits_ );
128  if ( fabs(corr_d0) < cut(indexD0_, double()) || ignoreCut(indexD0_) ) passCut(ret, indexD0_ );
129  if ( isConv || ignoreCut(indexConvRej_) ) passCut(ret, indexConvRej_ );
130  if ( pfIso < cut(indexPFIso_, double()) || ignoreCut(indexPFIso_) ) passCut(ret, indexPFIso_ );
131  if ( mva > cut(indexMVA_, double()) || ignoreCut(indexMVA_) ) passCut(ret, indexMVA_ );
132  if ( electronIDboolean || ignoreCut(indexElectronId_) ) passCut(ret, indexElectronId_);
133  setIgnored(ret);
134  return (bool)ret;
135  }
136 
137  private: // member variables
138 
140 
148 
150 };
151 
152 #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:106
virtual double et() const GCC11_FINAL
transverse energy
dictionary parameters
Definition: Parameters.py:2
bool spring11Cuts(const pat::Electron &electron, pat::strbitset &ret)
std::string electronIDvalue_
bool exists(std::string const &parameterName) const
checks if a parameter exists
reco::GsfTrackRef gsfTrack() const
override the reco::GsfElectron::gsfTrack method, to access the internal storage of the supercluster ...
Definition: Electron.cc:131
pat::strbitset::index_type index_type
Definition: Selector.h:30
void setIgnored(pat::strbitset &ret)
set ignored bits
Definition: Selector.h:225
pat::strbitset retInternal_
internal ret if users don&#39;t care about return bits
Definition: Selector.h:288
float convDist() const
Definition: GsfElectron.h:498
float convDcot() const
Definition: GsfElectron.h:499
float mva() const
Definition: GsfElectron.h:567
pat::strbitset bits_
the bitset indexed by strings
Definition: Selector.h:287
void passCut(pat::strbitset &ret, std::string const &s)
Passing cuts.
Definition: Selector.h:177
bool ignoreCut(std::string const &s) const
ignore the cut at index &quot;s&quot;
Definition: Selector.h:160
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:47
PFElectronSelector(edm::ParameterSet const &parameters)
Functor that operates on &lt;T&gt;
Definition: Selector.h:25
float electronID(const std::string &name) const
Returns a specific electron ID associated to the pat::Electron given its name.
Definition: Electron.cc:315
float userIsolation(IsolationKeys key) const
Definition: Lepton.h:51
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
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:213
void setIgnoredCuts(std::vector< std::string > const &bitsToIgnore)
set the bits to ignore from a vector
Definition: Selector.h:168
double dB(IpType type=None) const
Impact parameter wrt primary vertex or beamspot.
Definition: Electron.cc:378
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:195