Go to the documentation of this file.00001 #include "RecoTauTag/RecoTau/interface/RecoTauDiscriminantFunctions.h"
00002 #include "DataFormats/Math/interface/deltaR.h"
00003 #include "DataFormats/Math/interface/angle.h"
00004 #include <algorithm>
00005 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00006 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00007 #include "DataFormats/TauReco/interface/RecoTauPiZero.h"
00008 #include <boost/foreach.hpp>
00009
00010 namespace reco { namespace tau { namespace disc {
00011
00012
00013 namespace {
00014
00015 const PFCandidate& removeRef(const PFCandidateRef& pfRef) {
00016 return *pfRef;
00017 }
00018
00019 const RecoTauPiZero& removeRef(const RecoTauPiZero& piZero) {
00020 return piZero;
00021 }
00022
00023
00024 template<typename Collection, typename Function>
00025 VDouble extract(const Collection& cands, Function func) {
00026
00027 VDouble output;
00028 output.reserve(cands.size());
00029 for(typename Collection::const_iterator cand = cands.begin();
00030 cand != cands.end(); ++cand) {
00031 output.push_back(func(removeRef(*cand)));
00032 }
00033 return output;
00034 }
00035
00036 class DeltaRToAxis {
00037 public:
00038 DeltaRToAxis(const reco::Candidate::LorentzVector& axis):axis_(axis){}
00039 double operator()(const Candidate& cand)
00040 {
00041 return deltaR(cand.p4(), axis_);
00042 }
00043 private:
00044 const reco::Candidate::LorentzVector& axis_;
00045 };
00046
00047 }
00048
00049 PFCandidateRef mainTrack(Tau tau) {
00050 if (tau.signalPFChargedHadrCands().size() == 3) {
00051 for (size_t itrk = 0; itrk < 3; ++itrk) {
00052 if (tau.signalPFChargedHadrCands()[itrk]->charge() * tau.charge() < 0)
00053 return tau.signalPFChargedHadrCands()[itrk];
00054 }
00055 }
00056 return tau.leadPFChargedHadrCand();
00057 }
00058
00059 PFCandidateRefVector notMainTrack(Tau tau)
00060 {
00061 const PFCandidateRef& mainTrackRef = mainTrack(tau);
00062 PFCandidateRefVector output;
00063 output.reserve(tau.signalPFChargedHadrCands().size() - 1);
00064 BOOST_FOREACH(const PFCandidateRef& ref, tau.signalPFChargedHadrCands()) {
00065 if (ref != mainTrackRef)
00066 output.push_back(ref);
00067 }
00068 return output;
00069 }
00070
00071 double Pt(Tau tau) { return tau.pt(); }
00072 double Eta(Tau tau) { return tau.eta(); }
00073 double AbsEta(Tau tau) { return std::abs(tau.eta()); }
00074 double Mass(Tau tau) { return tau.mass(); }
00075 double DecayMode(Tau tau) { return tau.decayMode(); }
00076 double InvariantMassOfSignal(Tau tau) { return tau.mass(); }
00077
00078
00079
00080
00081
00082 double JetPt(Tau tau) {
00083 return tau.jetRef()->pt();
00084 }
00085
00086 double JetEta(Tau tau) {
00087 return tau.jetRef()->eta();
00088 }
00089
00090 double AbsJetEta(Tau tau) {
00091 return std::abs(tau.jetRef()->eta());
00092 }
00093
00094 double JetWidth(Tau tau) {
00095 return std::sqrt(
00096 std::abs(tau.jetRef()->etaetaMoment()) +
00097 std::abs(tau.jetRef()->phiphiMoment()));
00098 }
00099
00100 double JetTauDR(Tau tau) {
00101 return reco::deltaR(tau.p4(), tau.jetRef()->p4());
00102 }
00103
00104 double SignalPtFraction(Tau tau) {
00105 return tau.pt()/tau.jetRef()->pt();
00106 }
00107
00108 double IsolationChargedPtFraction(Tau tau) {
00109 return tau.isolationPFChargedHadrCandsPtSum()/tau.jetRef()->pt();
00110 }
00111
00112 double IsolationECALPtFraction(Tau tau) {
00113 return tau.isolationPFGammaCandsEtSum()/tau.jetRef()->pt();
00114 }
00115
00116 double IsolationNeutralHadronPtFraction(Tau tau) {
00117 double sum = 0.0;
00118 BOOST_FOREACH(PFCandidateRef cand, tau.isolationPFNeutrHadrCands()) {
00119 sum += cand->pt();
00120 }
00121 return sum/tau.jetRef()->pt();
00122 }
00123
00124 double ScaledEtaJetCollimation(Tau tau) {
00125 return tau.jetRef()->pt()*sqrt(std::abs(
00126 tau.jetRef()->etaetaMoment()));
00127 }
00128
00129 double OpeningDeltaR(Tau tau) {
00130 double sumEt = 0;
00131 double weightedDeltaR = 0;
00132 BOOST_FOREACH(const reco::PFCandidateRef& cand, tau.signalPFCands()) {
00133 double candEt = cand->et();
00134 double candDeltaR = reco::deltaR(cand->p4(), tau.p4());
00135 sumEt += candEt;
00136 weightedDeltaR += candDeltaR*candEt;
00137 }
00138 return (sumEt > 0) ? weightedDeltaR/sumEt : 0.0;
00139 }
00140
00141 double OpeningAngle3D(Tau tau) {
00142 double sumE = 0;
00143 double weightedAngle = 0;
00144 BOOST_FOREACH(const reco::PFCandidateRef& cand, tau.signalPFCands()) {
00145 double candE = cand->energy();
00146 double candAngle = angle(cand->p4(), tau.p4());
00147 sumE += candE;
00148 weightedAngle += candAngle*candE;
00149 }
00150 return (sumE > 0) ? weightedAngle/sumE : 0.0;
00151 }
00152
00153 double ScaledOpeningDeltaR(Tau tau) {
00154 double max = 0.0;
00155 const PFCandidateRefVector& cands = tau.signalPFCands();
00156 for (size_t i = 0; i < cands.size()-1; ++i) {
00157 for (size_t j = i+1; j < cands.size(); ++j) {
00158 double deltaRVal = deltaR(cands[i]->p4(), cands[j]->p4());
00159 if (deltaRVal > max) {
00160 max = deltaRVal;
00161 }
00162 }
00163 }
00164
00165 if ( max < 0.05 )
00166 max = 0.05;
00167
00168 return max*tau.pt();;
00169 }
00170
00171 double ScaledPhiJetCollimation(Tau tau) {
00172 return tau.jetRef()->pt()*sqrt(std::abs(
00173 tau.jetRef()->phiphiMoment()));
00174 }
00175
00176 double IsolationChargedAveragePtFraction(Tau tau) {
00177 size_t nIsoCharged = tau.isolationPFChargedHadrCands().size();
00178 double averagePt = (nIsoCharged) ?
00179 tau.isolationPFChargedHadrCandsPtSum()/nIsoCharged : 0;
00180 return averagePt/tau.leadPFChargedHadrCand()->pt();
00181 }
00182
00183 double MainTrackPtFraction(Tau tau) {
00184 return mainTrack(tau)->pt()/tau.jetRef()->pt();
00185 }
00186
00187 VDouble Dalitz2(Tau tau) {
00188 PFCandidateRef theMainTrack = mainTrack(tau);
00189 PFCandidateRefVector otherSignalTracks = notMainTrack(tau);
00190 const std::vector<RecoTauPiZero> &pizeros = tau.signalPiZeroCandidates();
00191 VDouble output;
00192 output.reserve(otherSignalTracks.size() + pizeros.size());
00193
00194 BOOST_FOREACH(PFCandidateRef trk, otherSignalTracks) {
00195 reco::Candidate::LorentzVector p4 = theMainTrack->p4() + trk->p4();
00196 output.push_back(p4.mass());
00197 }
00198
00199 BOOST_FOREACH(const RecoTauPiZero &pizero, pizeros) {
00200 reco::Candidate::LorentzVector p4 = theMainTrack->p4() + pizero.p4();
00201 output.push_back(p4.mass());
00202 }
00203 return output;
00204 }
00205
00206 double IsolationChargedSumHard(Tau tau) {
00207 VDouble isocands = extract(tau.isolationPFChargedHadrCands(),
00208 std::mem_fun_ref(&PFCandidate::pt));
00209 double output = 0.0;
00210 BOOST_FOREACH(double pt, isocands) {
00211 if (pt > 1.0)
00212 output += pt;
00213 }
00214 return output;
00215 }
00216
00217 double IsolationChargedSumSoft(Tau tau) {
00218 VDouble isocands = extract(tau.isolationPFChargedHadrCands(),
00219 std::mem_fun_ref(&PFCandidate::pt));
00220 double output = 0.0;
00221 BOOST_FOREACH(double pt, isocands) {
00222 if (pt < 1.0)
00223 output += pt;
00224 }
00225 return output;
00226 }
00227
00228
00229 double IsolationChargedSumHardRelative(Tau tau) {
00230 return IsolationChargedSumHard(tau)/tau.jetRef()->pt();
00231 }
00232
00233 double IsolationChargedSumSoftRelative(Tau tau) {
00234 return IsolationChargedSumSoft(tau)/tau.jetRef()->pt();
00235 }
00236
00237 double IsolationECALSumHard(Tau tau) {
00238 VDouble isocands = extract(tau.isolationPFGammaCands(),
00239 std::mem_fun_ref(&PFCandidate::pt));
00240 double output = 0.0;
00241 BOOST_FOREACH(double pt, isocands) {
00242 if (pt > 1.5)
00243 output += pt;
00244 }
00245 return output;
00246 }
00247
00248 double IsolationECALSumSoft(Tau tau) {
00249 VDouble isocands = extract(tau.isolationPFGammaCands(),
00250 std::mem_fun_ref(&PFCandidate::pt));
00251 double output = 0.0;
00252 BOOST_FOREACH(double pt, isocands) {
00253 if (pt < 1.5)
00254 output += pt;
00255 }
00256 return output;
00257 }
00258
00259
00260 double IsolationECALSumHardRelative(Tau tau) {
00261 return IsolationECALSumHard(tau)/tau.jetRef()->pt();
00262 }
00263 double IsolationECALSumSoftRelative(Tau tau) {
00264 return IsolationECALSumSoft(tau)/tau.jetRef()->pt();
00265 }
00266
00267 double EMFraction(Tau tau) {
00268 return tau.emFraction();
00269 }
00270
00271 double ImpactParameterSignificance(Tau tau) {
00272 return std::abs(tau.leadPFChargedHadrCandsignedSipt());
00273 }
00274
00275 double OutlierN(Tau tau) {
00276 return tau.isolationPFChargedHadrCands().size() +
00277 tau.isolationPFGammaCands().size();
00278 }
00279
00280 double OutlierNCharged(Tau tau) {
00281 return tau.isolationPFChargedHadrCands().size();
00282 }
00283
00284 double MainTrackPt(Tau tau) {
00285 PFCandidateRef trk = mainTrack(tau);
00286 return (!trk) ? 0.0 : trk->pt();
00287 }
00288
00289 double MainTrackEta(Tau tau) {
00290 PFCandidateRef trk = mainTrack(tau);
00291 return (!trk) ? 0.0 : trk->eta();
00292 }
00293
00294 double MainTrackAngle(Tau tau) {
00295 PFCandidateRef trk = mainTrack(tau);
00296 return (!trk) ? 0.0 : deltaR(trk->p4(), tau.p4());
00297 }
00298
00299 double OutlierSumPt(Tau tau) {
00300 return tau.isolationPFChargedHadrCandsPtSum() +
00301 tau.isolationPFGammaCandsEtSum();
00302 }
00303
00304 double ChargedOutlierSumPt(Tau tau) {
00305 return tau.isolationPFChargedHadrCandsPtSum();
00306 }
00307
00308 double NeutralOutlierSumPt(Tau tau) {
00309 return tau.isolationPFGammaCandsEtSum();
00310 }
00311
00312
00313 VDouble TrackPt(Tau tau) {
00314 return extract(notMainTrack(tau), std::mem_fun_ref(&PFCandidate::pt));
00315 }
00316
00317 VDouble TrackEta(Tau tau) {
00318 return extract(notMainTrack(tau), std::mem_fun_ref(&PFCandidate::eta));
00319 }
00320
00321 VDouble TrackAngle(Tau tau) {
00322 return extract(notMainTrack(tau), DeltaRToAxis(tau.p4()));
00323 }
00324
00325
00326 VDouble PiZeroPt(Tau tau) {
00327 return extract(tau.signalPiZeroCandidates(), std::mem_fun_ref(&RecoTauPiZero::pt));
00328 }
00329
00330 VDouble PiZeroEta(Tau tau) {
00331 return extract(tau.signalPiZeroCandidates(), std::mem_fun_ref(&RecoTauPiZero::eta));
00332 }
00333
00334 VDouble PiZeroAngle(Tau tau) {
00335 return extract(tau.signalPiZeroCandidates(), DeltaRToAxis(tau.p4()));
00336 }
00337
00338
00339 VDouble OutlierPt(Tau tau) {
00340 return extract(tau.isolationPFCands(), std::mem_fun_ref(&PFCandidate::pt));
00341 }
00342
00343 VDouble OutlierAngle(Tau tau) {
00344 return extract(tau.isolationPFCands(), DeltaRToAxis(tau.p4()));
00345 }
00346
00347 VDouble ChargedOutlierPt(Tau tau) {
00348 return extract(tau.isolationPFChargedHadrCands(),
00349 std::mem_fun_ref(&PFCandidate::pt));
00350 }
00351
00352 VDouble ChargedOutlierAngle(Tau tau) {
00353 return extract(tau.isolationPFChargedHadrCands(), DeltaRToAxis(tau.p4()));
00354 }
00355
00356 VDouble NeutralOutlierPt(Tau tau) {
00357 return extract(tau.isolationPFGammaCands(),
00358 std::mem_fun_ref(&PFCandidate::pt));
00359 }
00360
00361 VDouble NeutralOutlierAngle(Tau tau) {
00362 return extract(tau.isolationPFGammaCands(), DeltaRToAxis(tau.p4()));
00363 }
00364
00365
00366 VDouble Dalitz(Tau tau) {
00367 return Dalitz2(tau);
00368 }
00369
00370
00371
00372 VDouble FilteredObjectPt(Tau tau) { return VDouble(); }
00373 VDouble GammaOccupancy(Tau tau) { return VDouble(); }
00374 VDouble GammaPt(Tau tau) { return VDouble(); }
00375 VDouble InvariantMassOfSignalWithFiltered(Tau tau) { return VDouble(); }
00376 VDouble InvariantMass(Tau tau) { return VDouble(); }
00377 VDouble OutlierMass(Tau tau) { return VDouble(); }
00378
00379 }}}
00380