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