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 223 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.

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

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

Referenced by filterCandRef().

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

Definition at line 391 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().

391  {
392  switch(cand.particleId()) {
393  case PFCandidate::gamma:
394  return filterGammaCand(cand);
395  case PFCandidate::h0:
396  return filterNeutralHadronCand(cand);
397  // We use the same qcuts for muons/electrons and charged hadrons.
398  case PFCandidate::h:
399  case PFCandidate::e:
400  case PFCandidate::mu:
401  // no cuts ATM (track cuts applied in filterCand)
402  return true;
403  // Return false if we dont' know how to deal with this particle type
404  default:
405  return false;
406  };
407  return false;
408 }
bool filterNeutralHadronCand(const reco::PFCandidate &cand) const
bool filterGammaCand(const reco::PFCandidate &cand) const
virtual ParticleType particleId() const
Definition: PFCandidate.h:368
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 381 of file RecoTauQualityCuts.cc.

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

Referenced by filterCandByType().

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

Definition at line 386 of file RecoTauQualityCuts.cc.

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

Referenced by filterCandByType().

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

Filter a single Track.

Definition at line 335 of file RecoTauQualityCuts.cc.

References filterTrack_().

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

336 {
337  return filterTrack_(track);
338 }
bool filterTrack_(const T &trackRef) const
bool reco::tau::RecoTauQualityCuts::filterTrack ( const reco::TrackRef track) const

Definition at line 340 of file RecoTauQualityCuts.cc.

References filterTrack_().

341 {
342  return filterTrack_(track);
343 }
bool filterTrack_(const T &trackRef) const
template<typename T >
bool reco::tau::RecoTauQualityCuts::filterTrack_ ( const T trackRef) const
private

Definition at line 346 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().

347 {
348  const Track *track = trackRef.get();
349  if(minTrackPt_ >= 0 && !(track->pt() > minTrackPt_)) return false;
350  if(maxTrackChi2_ >= 0 && !(track->normalizedChi2() <= maxTrackChi2_)) return false;
351  if(checkHitPattern_) {
352  const reco::HitPattern &hitPattern = track->hitPattern();
353  if(minTrackPixelHits_ > 0 && !(hitPattern.numberOfValidPixelHits() >= minTrackPixelHits_)) return false;
354  if(minTrackHits_ > 0 && !(hitPattern.numberOfValidHits() >= minTrackHits_)) return false;
355  }
356  if(checkPV_ && pv_.isNull()) {
357  edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
358  "RecoTauQualityCuts is invalid. - filterTrack";
359  return false;
360  }
361 
363  !(std::fabs(track->dxy(pv_->position())) <= maxTransverseImpactParameter_))
364  return false;
365  if(maxDeltaZ_ >= 0 && !(std::fabs(track->dz(pv_->position())) <= maxDeltaZ_)) return false;
366  if(maxDeltaZToLeadTrack_ >= 0) {
367  if ( leadTrack_.isNull()) {
368  edm::LogError("QCutsNoValidLeadTrack") << "Lead track Ref in " <<
369  "RecoTauQualityCuts is invalid. - filterTrack";
370  return false;
371  }
372 
373  if(!(std::fabs(track->dz(pv_->position()) - leadTrack_->dz(pv_->position())) <= maxDeltaZToLeadTrack_))
374  return false;
375  }
376  if(minTrackVertexWeight_ > -1.0 && !(pv_->trackWeight(convertRef(trackRef)) >= minTrackVertexWeight_)) return false;
377 
378  return true;
379 }
int numberOfValidHits() const
Definition: HitPattern.h:737
uint16_t hitPattern[ARRAY_LENGTH]
Definition: HitPattern.h:410
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:585
bool isNull() const
Checks for null.
Definition: RefToBase.h:271
int numberOfValidPixelHits() const
Definition: HitPattern.h:752
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
JetCorrectorParametersCollection coll
Definition: classes.h:10
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 428 of file RecoTauQualityCuts.cc.

References leadTrack_.

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

Definition at line 433 of file RecoTauQualityCuts.cc.

References leadTrack_.

434 {
435  leadTrack_ = getTrackRef(leadCand);
436 }
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 438 of file RecoTauQualityCuts.cc.

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

439 {
440  if ( leadCand.isNonnull() ) {
441  leadTrack_ = getTrackRef(*leadCand);
442  } else {
443  // Set null
445  }
446 }
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:31
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().