CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

reco::tau::RecoTauQualityCuts Class Reference

#include <RecoTauQualityCuts.h>

List of all members.

Public Types

typedef boost::function< bool(const
PFCandidate &)> 
QCutFunc
typedef std::vector< QCutFuncQCutFuncCollection
typedef std::map
< PFCandidate::ParticleType,
QCutFuncCollection
QCutFuncMap

Public Member Functions

bool filter (const reco::PFCandidate &cand) const
 Filter a single PFCandidate.
template<typename PFCandRefType >
bool filterRef (const PFCandRefType &cand) const
 Filter a PFCandidate held by a smart pointer or Ref.
template<typename Coll >
Coll filterRefs (const Coll &refcoll, bool invert=false) const
 Filter a ref vector of PFCandidates.
const QCutFuncpredicate () const
 Get the predicate used to filter.
 RecoTauQualityCuts (const edm::ParameterSet &qcuts)
void setLeadTrack (const reco::PFCandidateRef &leadCand) const
void setLeadTrack (const reco::PFCandidate &leadCand) const
 Update the leading track.
void setPV (const reco::VertexRef &vtx) const
 Update the primary vertex.

Private Attributes

reco::TrackBaseRef leadTrack_
QCutFunc predicate_
reco::VertexRef pv_
QCutFuncMap qcuts_

Detailed Description

Definition at line 33 of file RecoTauQualityCuts.h.


Member Typedef Documentation

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

Definition at line 36 of file RecoTauQualityCuts.h.

Definition at line 37 of file RecoTauQualityCuts.h.

Definition at line 38 of file RecoTauQualityCuts.h.


Constructor & Destructor Documentation

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

Definition at line 142 of file RecoTauQualityCuts.cc.

References alignCSCRings::e, reco::tau::qcuts::etMin(), Exception, edm::ParameterSet::exists(), reco::PFCandidate::gamma, edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNames(), h, reco::PFCandidate::h0, leadTrack_, reco::tau::qcuts::mapAndCutByType(), reco::tau::qcuts::minTrackVertexWeight(), RPCpg::mu, predicate_, reco::tau::qcuts::ptMin(), pv_, qcuts_, AlCaHLTBitMon_QueryRunRegistry::string, reco::tau::qcuts::trkChi2(), reco::tau::qcuts::trkLongitudinalImpactParameter(), reco::tau::qcuts::trkLongitudinalImpactParameterWrtTrack(), reco::tau::qcuts::trkPixelHits(), reco::tau::qcuts::trkTrackerHits(), and reco::tau::qcuts::trkTransverseImpactParameter().

                                                                       {
  // Setup all of our predicates
  QCutFuncCollection chargedHadronCuts;
  QCutFuncCollection gammaCuts;
  QCutFuncCollection neutralHadronCuts;

  // Make sure there are no extra passed options
  std::set<std::string> passedOptionSet;
  std::vector<std::string> passedOptions = qcuts.getParameterNames();

  BOOST_FOREACH(const std::string& option, passedOptions) {
    passedOptionSet.insert(option);
  }

  // Build all the QCuts for tracks
  if (qcuts.exists("minTrackPt")) {
    chargedHadronCuts.push_back(
        boost::bind(qcuts::ptMin, _1,
                    qcuts.getParameter<double>("minTrackPt")));
    passedOptionSet.erase("minTrackPt");
  }

  if (qcuts.exists("maxTrackChi2")) {
    chargedHadronCuts.push_back(
        boost::bind(qcuts::trkChi2, _1,
                    qcuts.getParameter<double>("maxTrackChi2")));
    passedOptionSet.erase("maxTrackChi2");
  }

  if (qcuts.exists("minTrackPixelHits")) {
    chargedHadronCuts.push_back(boost::bind(
            qcuts::trkPixelHits, _1,
            qcuts.getParameter<uint32_t>("minTrackPixelHits")));
    passedOptionSet.erase("minTrackPixelHits");
  }

  if (qcuts.exists("minTrackHits")) {
    chargedHadronCuts.push_back(boost::bind(
            qcuts::trkTrackerHits, _1,
            qcuts.getParameter<uint32_t>("minTrackHits")));
    passedOptionSet.erase("minTrackHits");
  }

  // The impact parameter functions are bound to our member PV, since they
  // need it to compute the discriminant value.
  if (qcuts.exists("maxTransverseImpactParameter")) {
    chargedHadronCuts.push_back(boost::bind(
            qcuts::trkTransverseImpactParameter, _1, &pv_,
            qcuts.getParameter<double>("maxTransverseImpactParameter")));
    passedOptionSet.erase("maxTransverseImpactParameter");
  }

  if (qcuts.exists("maxDeltaZ")) {
    chargedHadronCuts.push_back(boost::bind(
            qcuts::trkLongitudinalImpactParameter, _1, &pv_,
            qcuts.getParameter<double>("maxDeltaZ")));
    passedOptionSet.erase("maxDeltaZ");
  }

  if (qcuts.exists("maxDeltaZToLeadTrack")) {
    chargedHadronCuts.push_back(boost::bind(
            qcuts::trkLongitudinalImpactParameterWrtTrack, _1, &leadTrack_,
            qcuts.getParameter<double>("maxDeltaZToLeadTrack")));
    passedOptionSet.erase("maxDeltaZToLeadTrack");
  }

  // Require tracks to contribute a minimum weight to the associated vertex.
  if (qcuts.exists("minTrackVertexWeight")) {
    chargedHadronCuts.push_back(boost::bind(
          qcuts::minTrackVertexWeight, _1, &pv_,
          qcuts.getParameter<double>("minTrackVertexWeight")));
    passedOptionSet.erase("minTrackVertexWeight");
  }

  // Build the QCuts for gammas
  if (qcuts.exists("minGammaEt")) {
    gammaCuts.push_back(boost::bind(
            qcuts::etMin, _1, qcuts.getParameter<double>("minGammaEt")));
    passedOptionSet.erase("minGammaEt");
  }

  // Build QCuts for netural hadrons
  if (qcuts.exists("minNeutralHadronEt")) {
    neutralHadronCuts.push_back(boost::bind(
            qcuts::etMin, _1,
            qcuts.getParameter<double>("minNeutralHadronEt")));
    passedOptionSet.erase("minNeutralHadronEt");
  }

  // Check if there are any remaining unparsed QCuts
  if (passedOptionSet.size()) {
    std::string unParsedOptions;
    bool thereIsABadParameter = false;
    BOOST_FOREACH(const std::string& option, passedOptionSet) {
      // Workaround for HLT - TODO FIXME
      if (option == "useTracksInsteadOfPFHadrons") {
        // Crash if true - no one should have this option enabled.
        if (qcuts.getParameter<bool>("useTracksInsteadOfPFHadrons")) {
          throw cms::Exception("DontUseTracksInQcuts")
            << "The obsolete exception useTracksInsteadOfPFHadrons "
            << "is set to true in the quality cut config." << std::endl;
        }
        continue;
      }

      // If we get to this point, there is a real unknown parameter
      thereIsABadParameter = true;

      unParsedOptions += option;
      unParsedOptions += "\n";
    }
    if (thereIsABadParameter) {
      throw cms::Exception("BadQualityCutConfig")
        << " The PSet passed to the RecoTauQualityCuts class had"
        << " the following unrecognized options: " << std::endl
        << unParsedOptions;
    }
  }

  // Make sure there are at least some quality cuts
  size_t nCuts = chargedHadronCuts.size() + gammaCuts.size()
    + neutralHadronCuts.size();
  if (!nCuts) {
    throw cms::Exception("BadQualityCutConfig")
      << " No options were passed to the quality cut class!" << std::endl;
  }

  // Map our QCut collections to the particle Ids they are associated to.
  qcuts_[PFCandidate::h] = chargedHadronCuts;
  qcuts_[PFCandidate::gamma] = gammaCuts;
  qcuts_[PFCandidate::h0] = neutralHadronCuts;
  // We use the same qcuts for muons/electrons and charged hadrons.
  qcuts_[PFCandidate::e] = chargedHadronCuts;
  qcuts_[PFCandidate::mu] = chargedHadronCuts;

  // Build a final level predicate that works on any PFCand
  predicate_ = boost::bind(qcuts::mapAndCutByType, _1, boost::cref(qcuts_));
}

Member Function Documentation

bool reco::tau::RecoTauQualityCuts::filter ( const reco::PFCandidate cand) const [inline]

Filter a single PFCandidate.

Definition at line 56 of file RecoTauQualityCuts.h.

References predicate_.

Referenced by filterRef(), and reco::tau::RecoTauEnergyRecoveryPlugin::operator()().

                                                   {
      return predicate_(cand);
    }
template<typename PFCandRefType >
bool reco::tau::RecoTauQualityCuts::filterRef ( const PFCandRefType &  cand) const [inline]

Filter a PFCandidate held by a smart pointer or Ref.

Definition at line 62 of file RecoTauQualityCuts.h.

References filter().

Referenced by filterRefs().

{ return filter(*cand); }
template<typename Coll >
Coll reco::tau::RecoTauQualityCuts::filterRefs ( const Coll &  refcoll,
bool  invert = false 
) const [inline]
const QCutFunc& reco::tau::RecoTauQualityCuts::predicate ( ) const [inline]

Get the predicate used to filter.

Definition at line 53 of file RecoTauQualityCuts.h.

References predicate_.

{ return predicate_; }
void reco::tau::RecoTauQualityCuts::setLeadTrack ( const reco::PFCandidate leadCand) const

Update the leading track.

Definition at line 300 of file RecoTauQualityCuts.cc.

References leadTrack_.

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

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

                                              {
  if (leadCand.isNonnull()) {
    leadTrack_ = getTrackRef(*leadCand);
  } else {
    // Set null
    leadTrack_ = reco::TrackBaseRef();
  }
}
void reco::tau::RecoTauQualityCuts::setPV ( const reco::VertexRef vtx) const [inline]

Member Data Documentation

Definition at line 79 of file RecoTauQualityCuts.h.

Referenced by RecoTauQualityCuts(), and setLeadTrack().

Definition at line 83 of file RecoTauQualityCuts.h.

Referenced by filter(), predicate(), and RecoTauQualityCuts().

Definition at line 77 of file RecoTauQualityCuts.h.

Referenced by RecoTauQualityCuts(), and setPV().

Definition at line 81 of file RecoTauQualityCuts.h.

Referenced by RecoTauQualityCuts().