CMS 3D CMS Logo

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

#include <RecoTauQualityCuts.h>

Public Types

typedef boost::function< bool(const
PFCandidate &)> 
CandQCutFunc
 
typedef std::vector< CandQCutFuncCandQCutFuncCollection
 
typedef std::map
< PFCandidate::ParticleType,
CandQCutFuncCollection
CandQCutFuncMap
 
typedef boost::function< bool(const
TrackBaseRef &)> 
TrackQCutFunc
 
typedef std::vector
< TrackQCutFunc
TrackQCutFuncCollection
 

Public Member Functions

bool filterCand (const reco::PFCandidate &cand) const
 Filter a single PFCandidate. More...
 
template<typename PFCandRefType >
bool filterCandRef (const PFCandRefType &cand) const
 Filter a PFCandidate held by a smart pointer or Ref. More...
 
template<typename Coll >
Coll filterCandRefs (const Coll &refcoll, bool invert=false) const
 Filter a ref vector of PFCandidates. More...
 
bool filterTrack (const reco::TrackBaseRef &track) const
 Filter a single Track. More...
 
bool filterTrack (const reco::TrackRef &track) const
 
template<typename Coll >
Coll filterTracks (const Coll &coll, bool invert=false) const
 Filter a collection of Tracks. More...
 
 RecoTauQualityCuts (const edm::ParameterSet &qcuts)
 
void setLeadTrack (const reco::TrackRef &leadTrack) const
 Update the leading track. More...
 
void setLeadTrack (const reco::PFCandidate &leadCand) const
 
void setLeadTrack (const reco::PFCandidateRef &leadCand) const
 
void setPV (const reco::VertexRef &vtx) const
 Update the primary vertex. More...
 

Private Member Functions

bool filterCandByType (const reco::PFCandidate &cand) const
 
bool filterGammaCand (const reco::PFCandidate &cand) const
 
bool filterNeutralHadronCand (const reco::PFCandidate &cand) const
 
template<typename T >
bool filterTrack_ (const T &trackRef) const
 

Private Attributes

bool checkHitPattern_
 
bool checkPV_
 
reco::TrackBaseRef leadTrack_
 
double maxDeltaZ_
 
double maxDeltaZToLeadTrack_
 
double maxTrackChi2_
 
double maxTransverseImpactParameter_
 
double minGammaEt_
 
double minNeutralHadronEt_
 
int minTrackHits_
 
int minTrackPixelHits_
 
double minTrackPt_
 
double minTrackVertexWeight_
 
reco::VertexRef pv_
 

Detailed Description

Definition at line 34 of file RecoTauQualityCuts.h.

Member Typedef Documentation

typedef boost::function<bool (const PFCandidate&)> reco::tau::RecoTauQualityCuts::CandQCutFunc

Definition at line 40 of file RecoTauQualityCuts.h.

Definition at line 41 of file RecoTauQualityCuts.h.

Definition at line 42 of file RecoTauQualityCuts.h.

typedef boost::function<bool (const TrackBaseRef&)> reco::tau::RecoTauQualityCuts::TrackQCutFunc

Definition at line 38 of file RecoTauQualityCuts.h.

Definition at line 39 of file RecoTauQualityCuts.h.

Constructor & Destructor Documentation

reco::tau::RecoTauQualityCuts::RecoTauQualityCuts ( const edm::ParameterSet qcuts)
explicit

Definition at line 226 of file RecoTauQualityCuts.cc.

References checkHitPattern_, checkPV_, edm::hlt::Exception, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNames(), maxDeltaZ_, maxDeltaZToLeadTrack_, maxTrackChi2_, maxTransverseImpactParameter_, minGammaEt_, minNeutralHadronEt_, minTrackHits_, minTrackPixelHits_, minTrackPt_, minTrackVertexWeight_, mergeVDriftHistosByStation::name, and AlCaHLTBitMon_QueryRunRegistry::string.

227 {
228  // Setup all of our predicates
229  CandQCutFuncCollection chargedHadronCuts;
230  CandQCutFuncCollection gammaCuts;
231  CandQCutFuncCollection neutralHadronCuts;
232 
233  // Make sure there are no extra passed options
234  std::set<std::string> passedOptionSet;
235  std::vector<std::string> passedOptions = qcuts.getParameterNames();
236 
237  BOOST_FOREACH(const std::string& option, passedOptions) {
238  passedOptionSet.insert(option);
239  }
240 
241  unsigned int nCuts = 0;
242  auto getDouble = [&qcuts, &passedOptionSet, &nCuts](const std::string& name) {
243  if(qcuts.exists(name)) {
244  ++nCuts;
245  passedOptionSet.erase(name);
246  return qcuts.getParameter<double>(name);
247  }
248  return -1.0;
249  };
250  auto getUint = [&qcuts, &passedOptionSet, &nCuts](const std::string& name) -> unsigned int {
251  if(qcuts.exists(name)) {
252  ++nCuts;
253  passedOptionSet.erase(name);
254  return qcuts.getParameter<unsigned int>(name);
255  }
256  return 0;
257  };
258 
259  // Build all the QCuts for tracks
260  minTrackPt_ = getDouble("minTrackPt");
261  maxTrackChi2_ = getDouble("maxTrackChi2");
262  minTrackPixelHits_ = getUint("minTrackPixelHits");
263  minTrackHits_ = getUint("minTrackHits");
264  maxTransverseImpactParameter_ = getDouble("maxTransverseImpactParameter");
265  maxDeltaZ_ = getDouble("maxDeltaZ");
266  maxDeltaZToLeadTrack_ = getDouble("maxDeltaZToLeadTrack");
267  // Require tracks to contribute a minimum weight to the associated vertex.
268  minTrackVertexWeight_ = getDouble("minTrackVertexWeight");
269 
270  // Use bit-wise & to avoid conditional code
273  (maxDeltaZ_ >= 0) ||
274  (maxDeltaZToLeadTrack_ >= 0) ||
275  (minTrackVertexWeight_ >= 0);
276 
277  // Build the QCuts for gammas
278  minGammaEt_ = getDouble("minGammaEt");
279 
280  // Build QCuts for netural hadrons
281  minNeutralHadronEt_ = getDouble("minNeutralHadronEt");
282 
283  // Check if there are any remaining unparsed QCuts
284  if ( passedOptionSet.size() ) {
285  std::string unParsedOptions;
286  bool thereIsABadParameter = false;
287  BOOST_FOREACH( const std::string& option, passedOptionSet ) {
288  // Workaround for HLT - TODO FIXME
289  if ( option == "useTracksInsteadOfPFHadrons" ) {
290  // Crash if true - no one should have this option enabled.
291  if ( qcuts.getParameter<bool>("useTracksInsteadOfPFHadrons") ) {
292  throw cms::Exception("DontUseTracksInQcuts")
293  << "The obsolete exception useTracksInsteadOfPFHadrons "
294  << "is set to true in the quality cut config." << std::endl;
295  }
296  continue;
297  }
298 
299  // If we get to this point, there is a real unknown parameter
300  thereIsABadParameter = true;
301 
302  unParsedOptions += option;
303  unParsedOptions += "\n";
304  }
305  if ( thereIsABadParameter ) {
306  throw cms::Exception("BadQualityCutConfig")
307  << " The PSet passed to the RecoTauQualityCuts class had"
308  << " the following unrecognized options: " << std::endl
309  << unParsedOptions;
310  }
311  }
312 
313  // Make sure there are at least some quality cuts
314  if ( !nCuts ) {
315  throw cms::Exception("BadQualityCutConfig")
316  << " No options were passed to the quality cut class!" << std::endl;
317  }
318 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< CandQCutFunc > CandQCutFuncCollection
std::vector< std::string > getParameterNames() const

Member Function Documentation

bool reco::tau::RecoTauQualityCuts::filterCand ( const reco::PFCandidate cand) const

Filter a single PFCandidate.

Definition at line 413 of file RecoTauQualityCuts.cc.

References filterCandByType(), filterTrack_(), reco::PFCandidate::gsfTrackRef(), query::result, and reco::PFCandidate::trackRef().

Referenced by filterCandRef().

414 {
415  auto trackRef = cand.trackRef();
416  bool result = true;
417  if(trackRef.isNonnull()) {
418  result = filterTrack_(trackRef);
419  }
420  else {
421  auto gsfTrackRef = cand.gsfTrackRef();
422  if(gsfTrackRef.isNonnull()) {
423  result = filterTrack_(gsfTrackRef);
424  }
425  }
426  if(result)
427  result = filterCandByType(cand);
428  return result;
429 }
bool filterCandByType(const reco::PFCandidate &cand) const
bool filterTrack_(const T &trackRef) const
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:429
tuple result
Definition: query.py:137
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:467
bool reco::tau::RecoTauQualityCuts::filterCandByType ( const reco::PFCandidate cand) const
private

Definition at line 394 of file RecoTauQualityCuts.cc.

References reco::PFCandidate::e, filterGammaCand(), filterNeutralHadronCand(), reco::PFCandidate::gamma, reco::PFCandidate::h, reco::PFCandidate::h0, reco::PFCandidate::mu, and reco::PFCandidate::particleId().

Referenced by filterCand().

394  {
395  switch(cand.particleId()) {
396  case PFCandidate::gamma:
397  return filterGammaCand(cand);
398  case PFCandidate::h0:
399  return filterNeutralHadronCand(cand);
400  // We use the same qcuts for muons/electrons and charged hadrons.
401  case PFCandidate::h:
402  case PFCandidate::e:
403  case PFCandidate::mu:
404  // no cuts ATM (track cuts applied in filterCand)
405  return true;
406  // Return false if we dont' know how to deal with this particle type
407  default:
408  return false;
409  };
410  return false;
411 }
bool filterNeutralHadronCand(const reco::PFCandidate &cand) const
bool filterGammaCand(const reco::PFCandidate &cand) const
virtual ParticleType particleId() const
Definition: PFCandidate.h:355
template<typename PFCandRefType >
bool reco::tau::RecoTauQualityCuts::filterCandRef ( const PFCandRefType &  cand) const
inline

Filter a PFCandidate held by a smart pointer or Ref.

Definition at line 77 of file RecoTauQualityCuts.h.

References filterCand().

Referenced by filterCandRefs().

77 { return filterCand(*cand); }
bool filterCand(const reco::PFCandidate &cand) const
Filter a single PFCandidate.
template<typename Coll >
Coll reco::tau::RecoTauQualityCuts::filterCandRefs ( const Coll &  refcoll,
bool  invert = false 
) const
inline
bool reco::tau::RecoTauQualityCuts::filterGammaCand ( const reco::PFCandidate cand) const
private

Definition at line 384 of file RecoTauQualityCuts.cc.

References reco::LeafCandidate::et(), and minGammaEt_.

Referenced by filterCandByType().

384  {
385  if(minGammaEt_ >= 0 && !(cand.et() > minGammaEt_)) return false;
386  return true;
387 }
virtual double et() const GCC11_FINAL
transverse energy
bool reco::tau::RecoTauQualityCuts::filterNeutralHadronCand ( const reco::PFCandidate cand) const
private

Definition at line 389 of file RecoTauQualityCuts.cc.

References reco::LeafCandidate::et(), and minNeutralHadronEt_.

Referenced by filterCandByType().

389  {
390  if(minNeutralHadronEt_ >= 0 && !(cand.et() > minNeutralHadronEt_)) return false;
391  return true;
392 }
virtual double et() const GCC11_FINAL
transverse energy
bool reco::tau::RecoTauQualityCuts::filterTrack ( const reco::TrackBaseRef track) const

Filter a single Track.

Definition at line 338 of file RecoTauQualityCuts.cc.

References filterTrack_().

Referenced by filterTracks(), and reco::tau::PFRecoTauChargedHadronFromTrackPlugin::operator()().

339 {
340  return filterTrack_(track);
341 }
bool filterTrack_(const T &trackRef) const
bool reco::tau::RecoTauQualityCuts::filterTrack ( const reco::TrackRef track) const

Definition at line 343 of file RecoTauQualityCuts.cc.

References filterTrack_().

344 {
345  return filterTrack_(track);
346 }
bool filterTrack_(const T &trackRef) const
template<typename T >
bool reco::tau::RecoTauQualityCuts::filterTrack_ ( const T trackRef) const
private

Definition at line 349 of file RecoTauQualityCuts.cc.

References checkHitPattern_, checkPV_, reco::TrackBase::dxy(), reco::TrackBase::dz(), reco::TrackBase::hitPattern(), edm::RefToBase< T >::isNull(), edm::Ref< C, T, F >::isNull(), leadTrack_, maxDeltaZ_, maxDeltaZToLeadTrack_, maxTrackChi2_, maxTransverseImpactParameter_, minTrackHits_, minTrackPixelHits_, minTrackPt_, minTrackVertexWeight_, reco::TrackBase::normalizedChi2(), reco::HitPattern::numberOfValidHits(), reco::HitPattern::numberOfValidPixelHits(), reco::TrackBase::pt(), and pv_.

Referenced by filterCand(), and filterTrack().

350 {
351  const Track *track = trackRef.get();
352  if(minTrackPt_ >= 0 && !(track->pt() > minTrackPt_)) return false;
353  if(maxTrackChi2_ >= 0 && !(track->normalizedChi2() <= maxTrackChi2_)) return false;
354  if(checkHitPattern_) {
355  const reco::HitPattern hitPattern = track->hitPattern();
356  if(minTrackPixelHits_ > 0 && !(hitPattern.numberOfValidPixelHits() >= minTrackPixelHits_)) return false;
357  if(minTrackHits_ > 0 && !(hitPattern.numberOfValidHits() >= minTrackHits_)) return false;
358  }
359  if(checkPV_ && pv_.isNull()) {
360  edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
361  "RecoTauQualityCuts is invalid. - filterTrack";
362  return false;
363  }
364 
366  !(std::fabs(track->dxy(pv_->position())) <= maxTransverseImpactParameter_))
367  return false;
368  if(maxDeltaZ_ >= 0 && !(std::fabs(track->dz(pv_->position())) <= maxDeltaZ_)) return false;
369  if(maxDeltaZToLeadTrack_ >= 0) {
370  if ( leadTrack_.isNull()) {
371  edm::LogError("QCutsNoValidLeadTrack") << "Lead track Ref in " <<
372  "RecoTauQualityCuts is invalid. - filterTrack";
373  return false;
374  }
375 
376  if(!(std::fabs(track->dz(pv_->position()) - leadTrack_->dz(pv_->position())) <= maxDeltaZToLeadTrack_))
377  return false;
378  }
379  if(minTrackVertexWeight_ > -1.0 && !(pv_->trackWeight(convertRef(trackRef)) >= minTrackVertexWeight_)) return false;
380 
381  return true;
382 }
int numberOfValidHits() const
Definition: HitPattern.h:589
bool isNull() const
Checks for null.
Definition: RefToBase.h:270
bool isNull() const
Checks for null.
Definition: Ref.h:247
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:125
int numberOfValidPixelHits() const
Definition: HitPattern.h:601
template<typename Coll >
Coll reco::tau::RecoTauQualityCuts::filterTracks ( const Coll &  coll,
bool  invert = false 
) const
inline

Filter a collection of Tracks.

Definition at line 63 of file RecoTauQualityCuts.h.

References filterTrack(), and convertSQLitetoXML_cfg::output.

64  {
65  Coll output;
66  BOOST_FOREACH( const typename Coll::value_type track, coll ) {
67  if ( filterTrack(track)^invert ) output.push_back(track);
68  }
69  return output;
70  }
Container::value_type value_type
bool filterTrack(const reco::TrackBaseRef &track) const
Filter a single Track.
void reco::tau::RecoTauQualityCuts::setLeadTrack ( const reco::TrackRef leadTrack) const

Update the leading track.

Definition at line 431 of file RecoTauQualityCuts.cc.

References leadTrack_.

432 {
433  leadTrack_ = reco::TrackBaseRef(leadTrack);
434 }
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:22
void reco::tau::RecoTauQualityCuts::setLeadTrack ( const reco::PFCandidate leadCand) const

Definition at line 436 of file RecoTauQualityCuts.cc.

References leadTrack_.

437 {
438  leadTrack_ = getTrackRef(leadCand);
439 }
void reco::tau::RecoTauQualityCuts::setLeadTrack ( const reco::PFCandidateRef leadCand) const

Update the leading track (using reference) If null, this will set the lead track ref null.

Definition at line 441 of file RecoTauQualityCuts.cc.

References edm::Ref< C, T, F >::isNonnull(), and leadTrack_.

442 {
443  if ( leadCand.isNonnull() ) {
444  leadTrack_ = getTrackRef(*leadCand);
445  } else {
446  // Set null
448  }
449 }
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:22
void reco::tau::RecoTauQualityCuts::setPV ( const reco::VertexRef vtx) const
inline

Member Data Documentation

bool reco::tau::RecoTauQualityCuts::checkHitPattern_
private

Definition at line 111 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

bool reco::tau::RecoTauQualityCuts::checkPV_
private

Definition at line 112 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

reco::TrackBaseRef reco::tau::RecoTauQualityCuts::leadTrack_
mutableprivate

Definition at line 99 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and setLeadTrack().

double reco::tau::RecoTauQualityCuts::maxDeltaZ_
private

Definition at line 106 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

double reco::tau::RecoTauQualityCuts::maxDeltaZToLeadTrack_
private

Definition at line 107 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

double reco::tau::RecoTauQualityCuts::maxTrackChi2_
private

Definition at line 102 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

double reco::tau::RecoTauQualityCuts::maxTransverseImpactParameter_
private

Definition at line 105 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

double reco::tau::RecoTauQualityCuts::minGammaEt_
private

Definition at line 109 of file RecoTauQualityCuts.h.

Referenced by filterGammaCand(), and RecoTauQualityCuts().

double reco::tau::RecoTauQualityCuts::minNeutralHadronEt_
private

Definition at line 110 of file RecoTauQualityCuts.h.

Referenced by filterNeutralHadronCand(), and RecoTauQualityCuts().

int reco::tau::RecoTauQualityCuts::minTrackHits_
private

Definition at line 104 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

int reco::tau::RecoTauQualityCuts::minTrackPixelHits_
private

Definition at line 103 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

double reco::tau::RecoTauQualityCuts::minTrackPt_
private

Definition at line 101 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

double reco::tau::RecoTauQualityCuts::minTrackVertexWeight_
private

Definition at line 108 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and RecoTauQualityCuts().

reco::VertexRef reco::tau::RecoTauQualityCuts::pv_
mutableprivate

Definition at line 97 of file RecoTauQualityCuts.h.

Referenced by filterTrack_(), and setPV().