CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/RecoTauTag/RecoTau/src/RecoTauDiscriminantFunctions.cc

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 // Helper functions
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 // Function = A PFTau member function
00024 template<typename Collection, typename Function>
00025 VDouble extract(const Collection& cands, Function func) {
00026   // #define CALL_MEMBER_FN(object,ptrToMember)  ((object).*(ptrToMember))
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 } // end helper functions
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  * HPStanc variables
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   // Correct for resolution
00165   if ( max < 0.05 )
00166     max = 0.05;
00167   // Make invariant of PT
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   // Add combos with tracks
00194   BOOST_FOREACH(PFCandidateRef trk, otherSignalTracks) {
00195     reco::Candidate::LorentzVector p4 = theMainTrack->p4() + trk->p4();
00196     output.push_back(p4.mass());
00197   }
00198   // Add combos with pizeros
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 // Relative versions.
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 // Relative versions.
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   //double result = tau.emFraction();
00269   reco::Candidate::LorentzVector gammaP4;
00270   BOOST_FOREACH(const reco::PFCandidateRef& gamma, tau.signalPFGammaCands()) {
00271     gammaP4 += gamma->p4();
00272   }
00273   double result = gammaP4.pt()/tau.pt();
00274 
00275   if (result > 0.99) {
00276     std::cout << "EM fraction = " << result
00277       << tau << std::endl;
00278     tau.dump(std::cout);
00279     std::cout << "charged" << std::endl;
00280     BOOST_FOREACH(const reco::PFCandidateRef cand, tau.signalPFChargedHadrCands()) {
00281       std::cout << " pt: " << cand->pt() << " type: " << cand->particleId() <<  " key: " << cand.key() << std::endl;
00282     }
00283     std::cout << "gammas" << std::endl;
00284     BOOST_FOREACH(const reco::PFCandidateRef cand, tau.signalPFGammaCands()) {
00285       std::cout << " pt: " << cand->pt() << " type: " << cand->particleId() <<  " key: " << cand.key() << std::endl;
00286     }
00287   }
00288   return result;
00289 }
00290 
00291 double ImpactParameterSignificance(Tau tau) {
00292   return std::abs(tau.leadPFChargedHadrCandsignedSipt());
00293 }
00294 
00295 double OutlierN(Tau tau) {
00296   return tau.isolationPFChargedHadrCands().size() +
00297       tau.isolationPFGammaCands().size();
00298 }
00299 
00300 double OutlierNCharged(Tau tau) {
00301   return tau.isolationPFChargedHadrCands().size();
00302 }
00303 
00304 double MainTrackPt(Tau tau) {
00305   PFCandidateRef trk = mainTrack(tau);
00306   return (!trk) ? 0.0 : trk->pt();
00307 }
00308 
00309 double MainTrackEta(Tau tau) {
00310   PFCandidateRef trk = mainTrack(tau);
00311   return (!trk) ? 0.0 : trk->eta();
00312 }
00313 
00314 double MainTrackAngle(Tau tau) {
00315   PFCandidateRef trk = mainTrack(tau);
00316   return (!trk) ? 0.0 : deltaR(trk->p4(), tau.p4());
00317 }
00318 
00319 double OutlierSumPt(Tau tau) {
00320   return tau.isolationPFChargedHadrCandsPtSum() +
00321       tau.isolationPFGammaCandsEtSum();
00322 }
00323 
00324 double ChargedOutlierSumPt(Tau tau) {
00325   return tau.isolationPFChargedHadrCandsPtSum();
00326 }
00327 
00328 double NeutralOutlierSumPt(Tau tau) {
00329   return tau.isolationPFGammaCandsEtSum();
00330 }
00331 
00332 // Quantities associated to tracks - that are not the main track
00333 VDouble TrackPt(Tau tau) {
00334   return extract(notMainTrack(tau), std::mem_fun_ref(&PFCandidate::pt));
00335 }
00336 
00337 VDouble TrackEta(Tau tau) {
00338   return extract(notMainTrack(tau), std::mem_fun_ref(&PFCandidate::eta));
00339 }
00340 
00341 VDouble TrackAngle(Tau tau) {
00342   return extract(notMainTrack(tau), DeltaRToAxis(tau.p4()));
00343 }
00344 
00345 // Quantities associated to PiZeros
00346 VDouble PiZeroPt(Tau tau) {
00347   return extract(tau.signalPiZeroCandidates(), std::mem_fun_ref(&RecoTauPiZero::pt));
00348 }
00349 
00350 VDouble PiZeroEta(Tau tau) {
00351   return extract(tau.signalPiZeroCandidates(), std::mem_fun_ref(&RecoTauPiZero::eta));
00352 }
00353 
00354 VDouble PiZeroAngle(Tau tau) {
00355   return extract(tau.signalPiZeroCandidates(), DeltaRToAxis(tau.p4()));
00356 }
00357 
00358 // Isolation quantities
00359 VDouble OutlierPt(Tau tau) {
00360   return extract(tau.isolationPFCands(), std::mem_fun_ref(&PFCandidate::pt));
00361 }
00362 
00363 VDouble OutlierAngle(Tau tau) {
00364   return extract(tau.isolationPFCands(), DeltaRToAxis(tau.p4()));
00365 }
00366 
00367 VDouble ChargedOutlierPt(Tau tau) {
00368   return extract(tau.isolationPFChargedHadrCands(),
00369                  std::mem_fun_ref(&PFCandidate::pt));
00370 }
00371 
00372 VDouble ChargedOutlierAngle(Tau tau) {
00373   return extract(tau.isolationPFChargedHadrCands(), DeltaRToAxis(tau.p4()));
00374 }
00375 
00376 VDouble NeutralOutlierPt(Tau tau) {
00377   return extract(tau.isolationPFGammaCands(),
00378                  std::mem_fun_ref(&PFCandidate::pt));
00379 }
00380 
00381 VDouble NeutralOutlierAngle(Tau tau) {
00382   return extract(tau.isolationPFGammaCands(), DeltaRToAxis(tau.p4()));
00383 }
00384 
00385 // Invariant mass of main track with other combinations
00386 VDouble Dalitz(Tau tau) {
00387   return Dalitz2(tau);
00388 }
00389 
00390 // The below functions are deprecated.
00391 // Not used, for backwards compatability
00392 VDouble FilteredObjectPt(Tau tau) { return VDouble(); }
00393 VDouble GammaOccupancy(Tau tau) { return VDouble(); }
00394 VDouble GammaPt(Tau tau) { return VDouble(); }
00395 VDouble InvariantMassOfSignalWithFiltered(Tau tau) { return VDouble(); }
00396 VDouble InvariantMass(Tau tau) { return VDouble(); }
00397 VDouble OutlierMass(Tau tau) { return VDouble(); }
00398 
00399 }}} // end reco::tau::disc namespace
00400