CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoTauTag/RecoTau/src/RecoTauIsolationMasking.cc

Go to the documentation of this file.
00001 #include <boost/foreach.hpp>
00002 #include "RecoTauTag/RecoTau/interface/RecoTauIsolationMasking.h"
00003 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00004 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00005 #include "DataFormats/Math/interface/deltaR.h"
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "RecoParticleFlow/PFClusterTools/interface/PFEnergyResolution.h"
00008 
00009 namespace reco { namespace tau {
00010 
00011 namespace {
00012 class DRSorter {
00013   public:
00014     DRSorter(const reco::Candidate::LorentzVector& axis):axis_(axis){}
00015     template<typename T>
00016     bool operator()(const T& t1, const T& t2) const {
00017       return reco::deltaR(t1->p4(), axis_) < reco::deltaR(t2->p4(), axis_);
00018     }
00019   private:
00020     const reco::Candidate::LorentzVector axis_;
00021 };
00022 // Sort by descending pt
00023 class PtSorter {
00024   public:
00025     PtSorter(){}
00026     template<typename T>
00027     bool operator()(const T& t1, const T& t2) const {
00028       return t1->pt() > t2->pt();
00029     }
00030 };
00031 
00032 // Check if an object is within DR of a track collection
00033 class MultiTrackDRFilter {
00034   public:
00035     MultiTrackDRFilter(double deltaR, const reco::PFCandidateRefVector& trks)
00036       :deltaR_(deltaR),tracks_(trks){}
00037     template <typename T>
00038     bool operator()(const T& t) const {
00039       BOOST_FOREACH(const reco::PFCandidateRef& trk, tracks_) {
00040         if (reco::deltaR(trk->p4(), t->p4()) < deltaR_)
00041           return true;
00042       }
00043       return false;
00044     }
00045   private:
00046     double deltaR_;
00047     const reco::PFCandidateRefVector& tracks_;
00048 };
00049 
00050 double square(double x) { return x*x; }
00051 
00052 template<typename T>
00053 reco::PFCandidateRefVector convertRefCollection(const T& coll) {
00054   reco::PFCandidateRefVector output;
00055   output.reserve(coll.size());
00056   BOOST_FOREACH(const reco::PFCandidateRef cand, coll) {
00057     output.push_back(cand);
00058   }
00059   return output;
00060 }
00061 }
00062 
00063 RecoTauIsolationMasking::RecoTauIsolationMasking(
00064     const edm::ParameterSet& pset):resolutions_(new PFEnergyResolution) {
00065   ecalCone_ = pset.getParameter<double>("ecalCone");
00066   hcalCone_ = pset.getParameter<double>("hcalCone");
00067   finalHcalCone_ = pset.getParameter<double>("finalHcalCone");
00068   maxSigmas_ = pset.getParameter<double>("maxSigmas");
00069 }
00070 
00071 // Need to explicitly define this in the .cc so we can use auto_ptr + forward
00072 RecoTauIsolationMasking::~RecoTauIsolationMasking(){}
00073 
00074 RecoTauIsolationMasking::IsoMaskResult
00075 RecoTauIsolationMasking::mask(const reco::PFTau& tau) const {
00076   IsoMaskResult output;
00077 
00078   typedef std::list<reco::PFCandidateRef> PFCandList;
00079   // Copy original iso collections.
00080   std::copy(tau.isolationPFGammaCands().begin(),
00081       tau.isolationPFGammaCands().end(), std::back_inserter(output.gammas));
00082   std::copy(tau.isolationPFNeutrHadrCands().begin(),
00083       tau.isolationPFNeutrHadrCands().end(),
00084       std::back_inserter(output.h0s));
00085 
00086   std::vector<PFCandList*> courses;
00087   courses.push_back(&(output.h0s));
00088   courses.push_back(&(output.gammas));
00089   // Mask using each one of the tracks
00090   BOOST_FOREACH(const reco::PFCandidateRef& track,
00091       tau.signalPFChargedHadrCands()) {
00092     double trackerEnergy = track->energy();
00093     double linkedEcalEnergy = track->ecalEnergy();
00094     double linkedHcalEnergy = track->hcalEnergy();
00095     math::XYZPointF posAtCalo = track->positionAtECALEntrance();
00096     // Get the linked calo energies & their errors
00097     double linkedSumErrSquared = 0;
00098     linkedSumErrSquared += square(
00099         resolutions_->getEnergyResolutionEm(linkedEcalEnergy, posAtCalo.eta()));
00100     linkedSumErrSquared += square(
00101         resolutions_->getEnergyResolutionHad(
00102           linkedHcalEnergy, posAtCalo.eta(), posAtCalo.phi()));
00103 
00104     // energyDelta is the difference between associated Calo and tracker energy
00105     double energyDelta = linkedEcalEnergy + linkedHcalEnergy - trackerEnergy;
00106 
00107     // Sort the neutral hadrons by DR to the track
00108     //DRSorter sorter(track->p4());
00109     PtSorter sorter;
00110 
00111     BOOST_FOREACH(PFCandList* course, courses) {
00112       // Sort by deltaR to the current track
00113       course->sort(sorter);
00114       PFCandList::iterator toEatIter = course->begin();
00115       // While there are still candidates to eat in this course and they are
00116       // within the cone.
00117       while (toEatIter != course->end()) {
00118         const reco::PFCandidate& toEat = **toEatIter;
00119         double toEatEnergy = toEat.energy();
00120         double toEatErrorSq = square(resolution(toEat));
00121         // Check if we can absorb this candidate into the track.
00122         if (inCone(*track, **toEatIter) &&
00123             (energyDelta + toEatEnergy)/std::sqrt(
00124               linkedSumErrSquared + toEatErrorSq) < maxSigmas_ ) {
00125           energyDelta += toEatEnergy;
00126           linkedSumErrSquared += toEatErrorSq;
00127           toEatIter = course->erase(toEatIter);
00128         } else {
00129           // otherwise skip to the next one
00130           ++toEatIter;
00131         }
00132       }
00133     }
00134   }
00135   // Filter out any final HCAL objects with in cones about the tracks.
00136   // This removes upward fluctuating HCAL objects
00137   if (finalHcalCone_ > 0) {
00138     MultiTrackDRFilter hcalFinalFilter(finalHcalCone_,
00139         tau.signalPFChargedHadrCands());
00140     std::remove_if(output.h0s.begin(), output.h0s.end(), hcalFinalFilter);
00141   }
00142   return output;
00143 }
00144 
00145 double RecoTauIsolationMasking::resolution(
00146     const reco::PFCandidate& cand) const {
00147   if (cand.particleId() == reco::PFCandidate::h0) {
00148     // NB for HCAL it returns relative energy
00149     return cand.energy()*resolutions_->getEnergyResolutionHad(cand.energy(),
00150         cand.eta(), cand.phi());
00151   } else if (cand.particleId() == reco::PFCandidate::gamma) {
00152     return resolutions_->getEnergyResolutionEm(cand.energy(), cand.eta());
00153   } else if (cand.particleId() == reco::PFCandidate::e) {
00154     // FIXME what is the electron resolution??
00155     return 0.15;
00156   } else {
00157     edm::LogWarning("IsoMask::res (bad pf id)")
00158       << "Unknown PF ID: " << cand.particleId();
00159   }
00160   return -1;
00161 }
00162 
00163 bool RecoTauIsolationMasking::inCone(const reco::PFCandidate& track,
00164     const reco::PFCandidate& cand) const {
00165   double openingDR = reco::deltaR(track.positionAtECALEntrance(), cand.p4());
00166   if (cand.particleId() == reco::PFCandidate::h0) {
00167     return (openingDR < hcalCone_);
00168   } else if (cand.particleId() == reco::PFCandidate::gamma ||
00169       cand.particleId() == reco::PFCandidate::e) {
00170     return openingDR < ecalCone_;
00171   } else {
00172     edm::LogWarning("IsoMask::inCone (bad pf id)")
00173       << "Unknown PF ID: " << cand.particleId()
00174       << " " <<  reco::PFCandidate::e;
00175   }
00176   return -1;
00177 }
00178 
00179 }} // end namespace reco::tau