CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFJetIDSelectionFunctor.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PatUtils_interface_PFJetIDSelectionFunctor_h
2 #define PhysicsTools_PatUtils_interface_PFJetIDSelectionFunctor_h
3 
4 
24 
26 
27 #include <TMath.h>
28 class PFJetIDSelectionFunctor : public Selector<pat::Jet> {
29 
30  public: // interface
31 
34 
36 
38  {
39  std::string versionStr = params.getParameter<std::string>("version");
40  std::string qualityStr = params.getParameter<std::string>("quality");
41 
42  if ( versionStr == "FIRSTDATA" )
44  else
46 
47  if ( qualityStr == "LOOSE") quality_ = LOOSE;
48  else if ( qualityStr == "TIGHT") quality_ = TIGHT;
49  else quality_ = LOOSE;
50 
51  push_back("CHF" );
52  push_back("NHF" );
53  push_back("CEF" );
54  push_back("NEF" );
55  push_back("NCH" );
56  push_back("nConstituents");
57 
58 
59  // Set some default cuts for LOOSE, TIGHT
60  if ( quality_ == LOOSE ) {
61  set("CHF", 0.0);
62  set("NHF", 0.99);
63  set("CEF", 0.99);
64  set("NEF", 0.99);
65  set("NCH", 0);
66  set("nConstituents", 1);
67  } else if ( quality_ == TIGHT ) {
68  set("CHF", 0.0);
69  set("NHF", 0.9);
70  set("CEF", 0.99);
71  set("NEF", 0.9);
72  set("NCH", 0);
73  set("nConstituents", 1);
74  }
75 
76 
77  // Now check the configuration to see if the user changed anything
78  if ( params.exists("CHF") ) set("CHF", params.getParameter<double>("CHF") );
79  if ( params.exists("NHF") ) set("NHF", params.getParameter<double>("NHF") );
80  if ( params.exists("CEF") ) set("CEF", params.getParameter<double>("CEF") );
81  if ( params.exists("NEF") ) set("NEF", params.getParameter<double>("NEF") );
82  if ( params.exists("NCH") ) set("NCH", params.getParameter<int> ("NCH") );
83  if ( params.exists("nConstuents") ) set("nConstituents", params.getParameter<int> ("nConstituents") );
84 
85 
86 
87  indexNConstituents_ = index_type (&bits_, "nConstituents");
88  indexNEF_ = index_type (&bits_, "NEF");
89  indexNHF_ = index_type (&bits_, "NHF");
90  indexCEF_ = index_type (&bits_, "CEF");
91  indexCHF_ = index_type (&bits_, "CHF");
92  indexNCH_ = index_type (&bits_, "NCH");
93 
95 
96  }
97 
98 
100  Quality_t quality ) :
101  version_(version), quality_(quality)
102  {
103 
104  push_back("CHF" );
105  push_back("NHF" );
106  push_back("CEF" );
107  push_back("NEF" );
108  push_back("NCH" );
109  push_back("nConstituents");
110 
111 
112  // Set some default cuts for LOOSE, TIGHT
113  if ( quality_ == LOOSE ) {
114  set("CHF", 0.0);
115  set("NHF", 0.99);
116  set("CEF", 0.99);
117  set("NEF", 0.99);
118  set("NCH", 0);
119  set("nConstituents", 1);
120  } else if ( quality_ == TIGHT ) {
121  set("CHF", 0.0);
122  set("NHF", 0.9);
123  set("CEF", 0.99);
124  set("NEF", 0.9);
125  set("NCH", 0);
126  set("nConstituents", 1);
127  }
128 
129 
130  indexNConstituents_ = index_type (&bits_, "nConstituents");
131  indexNEF_ = index_type (&bits_, "NEF");
132  indexNHF_ = index_type (&bits_, "NHF");
133  indexCEF_ = index_type (&bits_, "CEF");
134  indexCHF_ = index_type (&bits_, "CHF");
135  indexNCH_ = index_type (&bits_, "NCH");
136 
138  }
139 
140 
141  //
142  // Accessor from PAT jets
143  //
145  {
146  if ( version_ == FIRSTDATA ) {
147  if ( jet.currentJECLevel() == "Uncorrected" )
148  return firstDataCuts( jet, ret );
149  else
150  return firstDataCuts( jet.correctedJet("Uncorrected"), ret );
151  }
152  else {
153  return false;
154  }
155  }
157 
158  //
159  // Accessor from *CORRECTED* 4-vector, EMF, and Jet ID.
160  // This can be used with reco quantities.
161  //
162  bool operator()( reco::PFJet const & jet,
163  pat::strbitset & ret )
164  {
165  if ( version_ == FIRSTDATA ) return firstDataCuts( jet, ret );
166  else {
167  return false;
168  }
169  }
170 
171  //
172  // cuts based on craft 08 analysis.
173  //
174  bool firstDataCuts( reco::Jet const & jet,
175  pat::strbitset & ret)
176  {
177 
178  ret.set(false);
179 
180  // cache some variables
181  double chf = 0.0;
182  double nhf = 0.0;
183  double cef = 0.0;
184  double nef = 0.0;
185  int nch = 0;
186  int nconstituents = 0;
187 
188  // Have to do this because pat::Jet inherits from reco::Jet but not reco::PFJet
189  reco::PFJet const * pfJet = dynamic_cast<reco::PFJet const *>(&jet);
190  pat::Jet const * patJet = dynamic_cast<pat::Jet const *>(&jet);
191  reco::BasicJet const * basicJet = dynamic_cast<reco::BasicJet const *>(&jet);
192 
193  if ( patJet != 0 ) {
194 
195  if ( patJet->isPFJet() ) {
196  chf = patJet->chargedHadronEnergyFraction();
197  nhf = ( patJet->neutralHadronEnergy() + patJet->HFHadronEnergy() ) / patJet->energy();
198  cef = patJet->chargedEmEnergyFraction();
199  nef = patJet->neutralEmEnergyFraction();
200  nch = patJet->chargedMultiplicity();
201  nconstituents = patJet->numberOfDaughters();
202  }
203  // Handle the special case where this is a composed jet for
204  // subjet analyses
205  else if ( patJet->isBasicJet() ) {
206  double e_chf = 0.0;
207  double e_nhf = 0.0;
208  double e_cef = 0.0;
209  double e_nef = 0.0;
210  nch = 0;
211  nconstituents = 0;
212 
213  for ( reco::Jet::const_iterator ibegin = patJet->begin(),
214  iend = patJet->end(), isub = ibegin;
215  isub != iend; ++isub ) {
216  reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
217  e_chf += pfsub->chargedHadronEnergy();
218  e_nhf += (pfsub->neutralHadronEnergy() + pfsub->HFHadronEnergy());
219  e_cef += pfsub->chargedEmEnergy();
220  e_nef += pfsub->neutralEmEnergy();
221  nch += pfsub->chargedMultiplicity();
222  nconstituents += pfsub->numberOfDaughters();
223  }
224  double e = patJet->energy();
225  if ( e > 0.000001 ) {
226  chf = e_chf / e;
227  nhf = e_nhf / e;
228  cef = e_cef / e;
229  nef = e_nef / e;
230  } else {
231  chf = nhf = cef = nef = 0.0;
232  }
233  }
234  } // end if pat jet
235  else if ( pfJet != 0 ) {
236  chf = pfJet->chargedHadronEnergyFraction();
237  nhf = ( pfJet->neutralHadronEnergy() + pfJet->HFHadronEnergy() ) / pfJet->energy();
238  cef = pfJet->chargedEmEnergyFraction();
239  nef = pfJet->neutralEmEnergyFraction();
240  nch = pfJet->chargedMultiplicity();
241  nconstituents = pfJet->numberOfDaughters();
242  } // end if PF jet
243  // Handle the special case where this is a composed jet for
244  // subjet analyses
245  else if ( basicJet != 0 ) {
246  double e_chf = 0.0;
247  double e_nhf = 0.0;
248  double e_cef = 0.0;
249  double e_nef = 0.0;
250  nch = 0;
251  nconstituents = 0;
252 
253  for ( reco::Jet::const_iterator ibegin = basicJet->begin(),
254  iend = patJet->end(), isub = ibegin;
255  isub != iend; ++isub ) {
256  reco::PFJet const * pfsub = dynamic_cast<reco::PFJet const *>( &*isub );
257  e_chf += pfsub->chargedHadronEnergy();
258  e_nhf += (pfsub->neutralHadronEnergy() + pfsub->HFHadronEnergy());
259  e_cef += pfsub->chargedEmEnergy();
260  e_nef += pfsub->neutralEmEnergy();
261  nch += pfsub->chargedMultiplicity();
262  nconstituents += pfsub->numberOfDaughters();
263  }
264  double e = basicJet->energy();
265  if ( e > 0.000001 ) {
266  chf = e_chf / e;
267  nhf = e_nhf / e;
268  cef = e_cef / e;
269  nef = e_nef / e;
270  }
271  } // end if basic jet
272 
273 
274  // Cuts for all |eta|:
275  if ( ignoreCut(indexNConstituents_) || nconstituents > cut(indexNConstituents_, int() ) ) passCut( ret, indexNConstituents_);
276  if ( ignoreCut(indexNEF_) || ( nef < cut(indexNEF_, double()) ) ) passCut( ret, indexNEF_);
277  if ( ignoreCut(indexNHF_) || ( nhf < cut(indexNHF_, double()) ) ) passCut( ret, indexNHF_);
278  // Cuts for |eta| < 2.4:
279  if ( ignoreCut(indexCEF_) || ( cef < cut(indexCEF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCEF_);
280  if ( ignoreCut(indexCHF_) || ( chf > cut(indexCHF_, double()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexCHF_);
281  if ( ignoreCut(indexNCH_) || ( nch > cut(indexNCH_, int()) || std::abs(jet.eta()) > 2.4 ) ) passCut( ret, indexNCH_);
282 
283  setIgnored( ret );
284  return (bool)ret;
285  }
286 
287  private: // member variables
288 
291 
298 
299 };
300 
301 #endif
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
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:135
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:322
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:93
Base class for all types of Jets.
Definition: Jet.h:21
Jet correctedJet(const std::string &level, const std::string &flavor="none", const std::string &set="") const
Definition: Jet.cc:254
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: Jet.h:557
float chargedEmEnergyFraction() const
chargedEmEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:326
bool exists(std::string const &parameterName) const
checks if a parameter exists
#define abs(x)
Definition: mlp_lapack.h:159
pat::strbitset::index_type index_type
Definition: Selector.h:30
void setIgnored(pat::strbitset &ret)
set ignored bits
Definition: Selector.h:225
bool firstDataCuts(reco::Jet const &jet, pat::strbitset &ret)
pat::strbitset retInternal_
internal ret if users don&#39;t care about return bits
Definition: Selector.h:288
float HFHadronEnergy() const
HFHadronEnergy.
Definition: PFJet.h:111
Jets made from CaloTowers.
Definition: BasicJet.h:21
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:148
Jets made from PFObjects.
Definition: PFJet.h:22
virtual double eta() const
momentum pseudorapidity
PFJetIDSelectionFunctor(Version_t version, Quality_t quality)
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:143
PFJetIDSelectionFunctor(edm::ParameterSet const &params)
float neutralEmEnergyFraction() const
neutralEmEnergyFraction
Definition: PFJet.h:145
virtual double energy() const
energy
pat::strbitset bits_
the bitset indexed by strings
Definition: Selector.h:287
virtual size_t numberOfDaughters() const
number of daughters
void passCut(pat::strbitset &ret, std::string const &s)
Passing cuts.
Definition: Selector.h:177
float HFHadronEnergy() const
HFHadronEnergy.
Definition: Jet.h:342
bool ignoreCut(std::string const &s) const
ignore the cut at index &quot;s&quot;
Definition: Selector.h:160
bool operator()(const pat::Jet &jet, pat::strbitset &ret)
This provides the interface for base classes to select objects.
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:47
Functor that operates on &lt;T&gt;
Definition: Selector.h:25
PF Jet selector for pat::Jets.
bool isPFJet() const
check to see if the jet is a reco::PFJet
Definition: Jet.h:222
float chargedEmEnergyFraction() const
chargedEmEnergyFraction
Definition: PFJet.h:137
virtual size_t numberOfDaughters() const
Definition: Jet.h:409
bool operator()(reco::PFJet const &jet, pat::strbitset &ret)
virtual const_iterator begin() const
first daughter const_iterator
strbitset & set(bool val=true)
set method of all bits
Definition: strbitset.h:144
std::string currentJECLevel() const
return the name of the current step of jet energy corrections
Definition: Jet.h:118
Analysis-level calorimeter jet class.
Definition: Jet.h:67
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:213
virtual const_iterator end() const
last daughter const_iterator
bool isBasicJet() const
check to see if the jet is no more than a reco::BasicJet
Definition: Jet.h:224
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:95
float neutralEmEnergyFraction() const
neutralEmEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:328
int chargedMultiplicity() const
chargedMultiplicity
Definition: Jet.h:585
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:91
int cut(index_type const &i, int val) const
Access the int cut values at index &quot;s&quot;.
Definition: Selector.h:195