00001 #include "RecoTauTag/RecoTau/interface/RecoTauQualityCuts.h"
00002 #include "DataFormats/VertexReco/interface/Vertex.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00005 #include "DataFormats/TrackReco/interface/Track.h"
00006
00007 #include <boost/bind.hpp>
00008
00009 namespace reco { namespace tau {
00010
00011 namespace {
00012
00013 reco::Track const * getTrack(const PFCandidate& cand) {
00014 reco::Track const * ret = cand.trackRef().get();
00015 if (!ret) ret = cand.gsfTrackRef().get();
00016 return ret;
00017 }
00018
00019 const reco::TrackBaseRef getTrackRef(const PFCandidate& cand) {
00020 if (cand.trackRef().isNonnull())
00021 return reco::TrackBaseRef(cand.trackRef());
00022 else if (cand.gsfTrackRef().isNonnull()) {
00023 return reco::TrackBaseRef(cand.gsfTrackRef());
00024 }
00025 return reco::TrackBaseRef();
00026 }
00027 }
00028
00029
00030 namespace qcuts {
00031
00032 bool ptMin(const PFCandidate& cand, double cut) {
00033 return cand.pt() > cut;
00034 }
00035
00036 bool etMin(const PFCandidate& cand, double cut) {
00037 return cand.et() > cut;
00038 }
00039
00040 bool trkPixelHits(const PFCandidate& cand, int cut) {
00041
00042 auto trk = getTrack(cand);
00043 if (!trk) return false;
00044 return trk->hitPattern().numberOfValidPixelHits() >= cut;
00045 }
00046
00047 bool trkTrackerHits(const PFCandidate& cand, int cut) {
00048 auto trk = getTrack(cand);
00049 if (!trk) return false;
00050 return trk->hitPattern().numberOfValidHits() >= cut;
00051 }
00052
00053 bool trkTransverseImpactParameter(const PFCandidate& cand,
00054 const reco::VertexRef* pv,
00055 double cut) {
00056 if (pv->isNull()) {
00057 edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
00058 "RecoTauQualityCuts is invalid. - trkTransverseImpactParameter";
00059 return false;
00060 }
00061 auto trk = getTrack(cand);
00062 if (!trk) return false;
00063 return std::abs(trk->dxy((*pv)->position())) <= cut;
00064 }
00065
00066 bool trkLongitudinalImpactParameter(const PFCandidate& cand,
00067 const reco::VertexRef* pv,
00068 double cut) {
00069 if (pv->isNull()) {
00070 edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
00071 "RecoTauQualityCuts is invalid. - trkLongitudinalImpactParameter";
00072 return false;
00073 }
00074 auto trk = getTrack(cand);
00075 if (!trk) return false;
00076 double difference = std::abs(trk->dz((*pv)->position()));
00077
00078
00079 return difference <= cut;
00080 }
00081
00083 bool trkLongitudinalImpactParameterWrtTrack(const PFCandidate& cand,
00084 const reco::TrackBaseRef* trk, double cut) {
00085 if (trk->isNull()) {
00086 edm::LogError("QCutsNoValidLeadTrack") << "Lead track Ref in " <<
00087 "RecoTauQualityCuts is invalid. - trkLongitudinalImpactParameterWrtTrack";
00088 return false;
00089 }
00090 auto candTrk = getTrack(cand);
00091 if (!candTrk) return false;
00092 double difference = std::abs((*trk)->vz() - candTrk->vz());
00093 return difference <= cut;
00094 }
00095
00096
00097 bool minTrackVertexWeight(const PFCandidate& cand, const reco::VertexRef* pv,
00098 double cut) {
00099 if (pv->isNull()) {
00100 edm::LogError("QCutsNoPrimaryVertex") << "Primary vertex Ref in " <<
00101 "RecoTauQualityCuts is invalid. - minTrackVertexWeight";
00102 return false;
00103 }
00104 auto trk = getTrackRef(cand);
00105 if (!trk) return false;
00106 double weight = (*pv)->trackWeight(trk);
00107 return weight >= cut;
00108 }
00109
00110 bool trkChi2(const PFCandidate& cand, double cut) {
00111 auto trk = getTrack(cand);
00112 if (!trk) return false;
00113 return trk->normalizedChi2() <= cut;
00114 }
00115
00116
00117 bool AND(const PFCandidate& cand,
00118 const RecoTauQualityCuts::QCutFuncCollection& cuts) {
00119 BOOST_FOREACH(const RecoTauQualityCuts::QCutFunc& func, cuts) {
00120 if (!func(cand))
00121 return false;
00122 }
00123 return true;
00124 }
00125
00126
00127 bool mapAndCutByType(const PFCandidate& cand,
00128 const RecoTauQualityCuts::QCutFuncMap& funcMap) {
00129
00130 RecoTauQualityCuts::QCutFuncMap::const_iterator cuts =
00131 funcMap.find(cand.particleId());
00132
00133 if (cuts == funcMap.end())
00134 return false;
00135 else
00136
00137 return AND(cand, cuts->second);
00138 }
00139
00140 }
00141
00142 RecoTauQualityCuts::RecoTauQualityCuts(const edm::ParameterSet &qcuts) {
00143
00144 QCutFuncCollection chargedHadronCuts;
00145 QCutFuncCollection gammaCuts;
00146 QCutFuncCollection neutralHadronCuts;
00147
00148
00149 std::set<std::string> passedOptionSet;
00150 std::vector<std::string> passedOptions = qcuts.getParameterNames();
00151
00152 BOOST_FOREACH(const std::string& option, passedOptions) {
00153 passedOptionSet.insert(option);
00154 }
00155
00156
00157 if (qcuts.exists("minTrackPt")) {
00158 chargedHadronCuts.push_back(
00159 boost::bind(qcuts::ptMin, _1,
00160 qcuts.getParameter<double>("minTrackPt")));
00161 passedOptionSet.erase("minTrackPt");
00162 }
00163
00164 if (qcuts.exists("maxTrackChi2")) {
00165 chargedHadronCuts.push_back(
00166 boost::bind(qcuts::trkChi2, _1,
00167 qcuts.getParameter<double>("maxTrackChi2")));
00168 passedOptionSet.erase("maxTrackChi2");
00169 }
00170
00171 if (qcuts.exists("minTrackPixelHits")) {
00172 chargedHadronCuts.push_back(boost::bind(
00173 qcuts::trkPixelHits, _1,
00174 qcuts.getParameter<uint32_t>("minTrackPixelHits")));
00175 passedOptionSet.erase("minTrackPixelHits");
00176 }
00177
00178 if (qcuts.exists("minTrackHits")) {
00179 chargedHadronCuts.push_back(boost::bind(
00180 qcuts::trkTrackerHits, _1,
00181 qcuts.getParameter<uint32_t>("minTrackHits")));
00182 passedOptionSet.erase("minTrackHits");
00183 }
00184
00185
00186
00187 if (qcuts.exists("maxTransverseImpactParameter")) {
00188 chargedHadronCuts.push_back(boost::bind(
00189 qcuts::trkTransverseImpactParameter, _1, &pv_,
00190 qcuts.getParameter<double>("maxTransverseImpactParameter")));
00191 passedOptionSet.erase("maxTransverseImpactParameter");
00192 }
00193
00194 if (qcuts.exists("maxDeltaZ")) {
00195 chargedHadronCuts.push_back(boost::bind(
00196 qcuts::trkLongitudinalImpactParameter, _1, &pv_,
00197 qcuts.getParameter<double>("maxDeltaZ")));
00198 passedOptionSet.erase("maxDeltaZ");
00199 }
00200
00201 if (qcuts.exists("maxDeltaZToLeadTrack")) {
00202 chargedHadronCuts.push_back(boost::bind(
00203 qcuts::trkLongitudinalImpactParameterWrtTrack, _1, &leadTrack_,
00204 qcuts.getParameter<double>("maxDeltaZToLeadTrack")));
00205 passedOptionSet.erase("maxDeltaZToLeadTrack");
00206 }
00207
00208
00209 if (qcuts.exists("minTrackVertexWeight")) {
00210 chargedHadronCuts.push_back(boost::bind(
00211 qcuts::minTrackVertexWeight, _1, &pv_,
00212 qcuts.getParameter<double>("minTrackVertexWeight")));
00213 passedOptionSet.erase("minTrackVertexWeight");
00214 }
00215
00216
00217 if (qcuts.exists("minGammaEt")) {
00218 gammaCuts.push_back(boost::bind(
00219 qcuts::etMin, _1, qcuts.getParameter<double>("minGammaEt")));
00220 passedOptionSet.erase("minGammaEt");
00221 }
00222
00223
00224 if (qcuts.exists("minNeutralHadronEt")) {
00225 neutralHadronCuts.push_back(boost::bind(
00226 qcuts::etMin, _1,
00227 qcuts.getParameter<double>("minNeutralHadronEt")));
00228 passedOptionSet.erase("minNeutralHadronEt");
00229 }
00230
00231
00232 if (passedOptionSet.size()) {
00233 std::string unParsedOptions;
00234 bool thereIsABadParameter = false;
00235 BOOST_FOREACH(const std::string& option, passedOptionSet) {
00236
00237 if (option == "useTracksInsteadOfPFHadrons") {
00238
00239 if (qcuts.getParameter<bool>("useTracksInsteadOfPFHadrons")) {
00240 throw cms::Exception("DontUseTracksInQcuts")
00241 << "The obsolete exception useTracksInsteadOfPFHadrons "
00242 << "is set to true in the quality cut config." << std::endl;
00243 }
00244 continue;
00245 }
00246
00247
00248 thereIsABadParameter = true;
00249
00250 unParsedOptions += option;
00251 unParsedOptions += "\n";
00252 }
00253 if (thereIsABadParameter) {
00254 throw cms::Exception("BadQualityCutConfig")
00255 << " The PSet passed to the RecoTauQualityCuts class had"
00256 << " the following unrecognized options: " << std::endl
00257 << unParsedOptions;
00258 }
00259 }
00260
00261
00262 size_t nCuts = chargedHadronCuts.size() + gammaCuts.size()
00263 + neutralHadronCuts.size();
00264 if (!nCuts) {
00265 throw cms::Exception("BadQualityCutConfig")
00266 << " No options were passed to the quality cut class!" << std::endl;
00267 }
00268
00269
00270 qcuts_[PFCandidate::h] = chargedHadronCuts;
00271 qcuts_[PFCandidate::gamma] = gammaCuts;
00272 qcuts_[PFCandidate::h0] = neutralHadronCuts;
00273
00274 qcuts_[PFCandidate::e] = chargedHadronCuts;
00275 qcuts_[PFCandidate::mu] = chargedHadronCuts;
00276
00277
00278 predicate_ = boost::bind(qcuts::mapAndCutByType, _1, boost::cref(qcuts_));
00279 }
00280
00281 std::pair<edm::ParameterSet, edm::ParameterSet> factorizePUQCuts(
00282 const edm::ParameterSet& input) {
00283
00284 edm::ParameterSet puCuts;
00285 edm::ParameterSet nonPUCuts;
00286
00287 std::vector<std::string> inputNames = input.getParameterNames();
00288 BOOST_FOREACH(const std::string& cut, inputNames) {
00289 if (cut == "minTrackVertexWeight" || cut == "maxDeltaZ"
00290 || cut == "maxDeltaZToLeadTrack") {
00291 puCuts.copyFrom(input, cut);
00292 } else {
00293 nonPUCuts.copyFrom(input, cut);
00294 }
00295 }
00296 return std::make_pair(puCuts, nonPUCuts);
00297 }
00298
00299
00300 void RecoTauQualityCuts::setLeadTrack(
00301 const reco::PFCandidate& leadCand) const {
00302 leadTrack_ = getTrackRef(leadCand);
00303 }
00304
00305 void RecoTauQualityCuts::setLeadTrack(
00306 const reco::PFCandidateRef& leadCand) const {
00307 if (leadCand.isNonnull()) {
00308 leadTrack_ = getTrackRef(*leadCand);
00309 } else {
00310
00311 leadTrack_ = reco::TrackBaseRef();
00312 }
00313 }
00314
00315
00316 }}