CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/RecoTauTag/RecoTau/src/RecoTauCommonUtilities.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00002 
00003 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00004 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00005 #include "DataFormats/JetReco/interface/PFJet.h"
00006 
00007 #include <algorithm>
00008 
00009 typedef std::vector<reco::PFCandidatePtr> PFCandPtrs;
00010 typedef PFCandPtrs::iterator PFCandIter;
00011 
00012 namespace reco { namespace tau {
00013 
00014 class SortPFCandsDescendingPt {
00015   public:
00016     bool operator()(const PFCandidatePtr& a, const PFCandidatePtr& b) const {
00017       return a->pt() > b->pt();
00018     }
00019 };
00020 
00021 std::vector<PFCandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>& piZeros) {
00022   std::vector<PFCandidatePtr> output;
00023 
00024   for(std::vector<RecoTauPiZero>::const_iterator piZero = piZeros.begin();
00025       piZero != piZeros.end(); ++piZero)
00026   {
00027     for(size_t iDaughter = 0; iDaughter < piZero->numberOfDaughters(); ++iDaughter)
00028     {
00029       output.push_back(PFCandidatePtr(piZero->daughterPtr(iDaughter)));
00030     }
00031   }
00032   return output;
00033 }
00034 
00035 std::vector<reco::PFCandidatePtr> pfCandidates(const reco::PFJet& jet, int particleId, bool sort)
00036 {
00037   PFCandPtrs pfCands = jet.getPFConstituents();
00038   PFCandPtrs selectedPFCands;
00039 
00040   for(PFCandIter cand = pfCands.begin(); cand != pfCands.end(); ++cand)
00041   {
00042     if( (**cand).particleId() == particleId )
00043       selectedPFCands.push_back(*cand);
00044   }
00045 
00046   if ( sort ) std::sort(selectedPFCands.begin(), selectedPFCands.end(), SortPFCandsDescendingPt());
00047 
00048   return selectedPFCands;
00049 }
00050 
00051 std::vector<reco::PFCandidatePtr> pfCandidates(const reco::PFJet& jet, const std::vector<int>& particleIds, bool sort)
00052 {
00053   PFCandPtrs output;
00054   // Get each desired candidate type, unsorted for now
00055   for(std::vector<int>::const_iterator particleId = particleIds.begin();
00056       particleId != particleIds.end(); ++particleId) {
00057     PFCandPtrs selectedPFCands = pfCandidates(jet, *particleId, false);
00058     output.insert(output.end(), selectedPFCands.begin(), selectedPFCands.end());
00059   }
00060   if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
00061   return output;
00062 }
00063 
00064 std::vector<reco::PFCandidatePtr> pfGammas(const reco::PFJet& jet, bool sort) {
00065   return pfCandidates(jet, reco::PFCandidate::gamma, sort);
00066 }
00067 
00068 std::vector<reco::PFCandidatePtr> pfChargedCands(const reco::PFJet& jet,
00069                                                  bool sort) {
00070   PFCandPtrs output;
00071   PFCandPtrs mus = pfCandidates(jet, reco::PFCandidate::mu, false);
00072   PFCandPtrs es = pfCandidates(jet, reco::PFCandidate::e, false);
00073   PFCandPtrs chs = pfCandidates(jet, reco::PFCandidate::h, false);
00074   output.reserve(mus.size() + es.size() + chs.size());
00075   output.insert(output.end(), mus.begin(), mus.end());
00076   output.insert(output.end(), es.begin(), es.end());
00077   output.insert(output.end(), chs.begin(), chs.end());
00078   if (sort) std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
00079   return output;
00080 }
00081 
00082 
00083 
00084 } }