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 
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  rhoTag_ = parameters.getParameter< edm::InputTag> ("rhoSrc");
29 
31 
32  if ( versionStr == "TOPPAG" ) {
33  version = TOPPAG;
34  }
35  else {
36  throw cms::Exception("InvalidInput") << "Expect version to be one of : TOPPAG" << std::endl;
37  }
38 
39 
40  initialize( version,
41  parameters.getParameter<bool> ("Fiducial"),
42  parameters.getParameter<int> ("MaxMissingHits"),
43  parameters.getParameter<double>("D0"),
44  parameters.getParameter<bool> ("ConversionRejection"),
45  parameters.getParameter<double>("PFIso"),
46  parameters.getParameter<double>("MVA")
47  );
48  if ( parameters.exists("cutsToIgnore") )
49  setIgnoredCuts( parameters.getParameter<std::vector<std::string> >("cutsToIgnore") );
50 
52 
53  }
54 
56  bool fid=true,
57  int nMissingHits = 0,
58  double d0 = 0.02,
59  bool convRej = true,
60  double pfiso = 0.10,
61  double mva = 0.4
62  )
63  {
64  version_ = version;
65 
66  // size_t found;
67  // found = eidUsed.find("NONE");
68  // if ( found != string::npos)
69 
70  push_back("Fiducial" );
71  push_back("MaxMissingHits", nMissingHits );
72  push_back("D0", d0 );
73  push_back("ConversionRejection" );
74  push_back("PFIso", pfiso );
75  push_back("MVA", mva );
76 
77  set("Fiducial", fid);
78  set("MaxMissingHits");
79  set("D0");
80  set("ConversionRejection", convRej);
81  set("PFIso");
82  set("MVA");
83 
84  indexFid_ = index_type(&bits_, "Fiducial" );
85  indexMaxMissingHits_ = index_type(&bits_, "MaxMissingHits" );
86  indexD0_ = index_type(&bits_, "D0" );
87  indexConvRej_ = index_type(&bits_, "ConversionRejection" );
88  indexPFIso_ = index_type(&bits_, "PFIso" );
89  indexMVA_ = index_type(&bits_, "MVA" );
90  }
91 
92  // Allow for multiple definitions of the cuts.
94  {
95  if (version_ == TOPPAG ) return topPAGRefCuts(electron, event, ret);
96  else {
97  return false;
98  }
99  }
100 
102 
103  // cuts based on top group L+J synchronization exercise
105  {
106 
107 
108 
109  ret.set(false);
110 
111  edm::Handle<double> rhoHandle;
112  event.getByLabel(rhoTag_, rhoHandle);
113 
114  double rhoIso = std::max(*(rhoHandle.product()), 0.0);
115 
116  double scEta = electron.superCluster()->eta();
117  double dB = electron.dB();
119  double chIso = electron.userIsolation(pat::PfChargedHadronIso);
120  double nhIso = electron.userIsolation(pat::PfNeutralHadronIso);
121  double phIso = electron.userIsolation(pat::PfGammaIso);
122  double pfIso = ( chIso + std::max(0.0, nhIso + phIso - rhoIso*AEff) )/ electron.ecalDrivenMomentum().pt();
123  int mHits = electron.gsfTrack()->trackerExpectedHitsInner().numberOfHits();
124  double mva = electron.electronID("mvaTrigV0");
125  //Electron Selection for e+jets
126  //-----------------------------
127  bool fid = ! (fabs(scEta) > 1.4442 && fabs(scEta) < 1.5660 );
128 
129 
130 
131  if ( fid || ignoreCut(indexFid_) ) passCut(ret, indexFid_ );
133  if ( fabs(dB) < cut(indexD0_, double()) || ignoreCut(indexD0_) ) passCut(ret, indexD0_ );
134  if ( electron.passConversionVeto() || ignoreCut(indexConvRej_) ) passCut(ret, indexConvRej_ );
135  if ( pfIso < cut(indexPFIso_, double()) || ignoreCut(indexPFIso_) ) passCut(ret, indexPFIso_ );
136  if ( mva > cut(indexMVA_, double()) || ignoreCut(indexMVA_) ) passCut(ret, indexMVA_ );
137  setIgnored(ret);
138  return (bool)ret;
139  }
140 
141  private: // member variables
142 
144 
151 
152 
153  std::string electronIDvalue_;
155 
156  bool operator()( const pat::Electron & electron, pat::strbitset & ret ) { return false;}
157 };
158 
159 #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
dictionary parameters
Definition: Parameters.py:2
const LorentzVector & ecalDrivenMomentum() const
Definition: Electron.h:187
reco::SuperClusterRef superCluster() const
override the reco::GsfElectron::superCluster method, to access the internal storage of the superclust...
Definition: Electron.cc:150
std::string electronIDvalue_
static Double_t GetElectronEffectiveArea(ElectronEffectiveAreaType type, Double_t SCEta, ElectronEffectiveAreaTarget EffectiveAreaTarget=kEleEAData2011)
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
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
const T & max(const T &a, const T &b)
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)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
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
T const * product() const
Definition: Handle.h:74
bool passConversionVeto() const
vertex fit combined with missing number of hits method
Definition: Electron.h:235
bool topPAGRefCuts(const pat::Electron &electron, edm::EventBase const &event, pat::strbitset &ret)
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:213
bool operator()(const pat::Electron &electron, edm::EventBase const &event, pat::strbitset &ret)
This provides an alternative signature that includes extra information.
void setIgnoredCuts(std::vector< std::string > const &bitsToIgnore)
set the bits to ignore from a vector
Definition: Selector.h:168
void initialize(Version_t version, bool fid=true, int nMissingHits=0, double d0=0.02, bool convRej=true, double pfiso=0.10, double mva=0.4)
list fid
Definition: NewTree.py:51
double dB(IpType type=None) const
Impact parameter wrt primary vertex or beamspot.
Definition: Electron.cc:378
int cut(index_type const &i, int val) const
Access the int cut values at index &quot;s&quot;.
Definition: Selector.h:195