CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTauTag/RecoTau/plugins/RecoTauPhotonFilter.cc

Go to the documentation of this file.
00001 /*
00002  * ===========================================================================
00003  *
00004  *       Filename:  RecoTauPhotonFilter
00005  *
00006  *    Description:  Modify taus to recursively filter lowpt photons
00007  *
00008  *         Author:  Evan K. Friis (UC Davis)
00009  *
00010  * ===========================================================================
00011  */
00012 
00013 #include <boost/foreach.hpp>
00014 #include "RecoTauTag/RecoTau/interface/RecoTauBuilderPlugins.h"
00015 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00016 
00017 namespace reco { namespace tau {
00018 
00019 // Filter photons
00020 class RecoTauPhotonFilter : public RecoTauModifierPlugin {
00021   public:
00022     explicit RecoTauPhotonFilter(const edm::ParameterSet& pset);
00023     virtual ~RecoTauPhotonFilter() {}
00024     void operator()(PFTau&) const;
00025   private:
00026     bool filter(const RecoTauPiZero* piZero,
00027                 const reco::Candidate::LorentzVector& tau) const;
00028     double minPtFractionSinglePhotons_;
00029     double minPtFractionPiZeroes_;
00030 };
00031 
00032 RecoTauPhotonFilter::RecoTauPhotonFilter(
00033     const edm::ParameterSet& pset):RecoTauModifierPlugin(pset) {
00034   minPtFractionSinglePhotons_ =
00035       pset.getParameter<double>("minPtFractionSinglePhotons");
00036   minPtFractionPiZeroes_ =
00037       pset.getParameter<double>("minPtFractionPiZeroes");
00038 }
00039 
00040 // Sort container of PiZero pointers by ascending Pt
00041 namespace {
00042 struct PiZeroPtSorter {
00043   bool operator()(const RecoTauPiZero* A, const RecoTauPiZero* B) {
00044     return (A->pt() < B->pt());
00045   }
00046 };
00047 }
00048 
00049 // Decide if we should filter a pi zero or not.
00050 bool RecoTauPhotonFilter::filter( const RecoTauPiZero* piZero,
00051     const reco::Candidate::LorentzVector& total) const {
00052   if (piZero->numberOfDaughters() > 1)
00053     return piZero->pt()/total.pt() < minPtFractionPiZeroes_;
00054   return piZero->pt()/total.pt() < minPtFractionSinglePhotons_;
00055 }
00056 
00057 void RecoTauPhotonFilter::operator()(PFTau& tau) const {
00058   std::vector<const RecoTauPiZero*> signalPiZeros;
00059   BOOST_FOREACH(const RecoTauPiZero& piZero, tau.signalPiZeroCandidates()) {
00060     signalPiZeros.push_back(&piZero);
00061   }
00062   std::sort(signalPiZeros.begin(), signalPiZeros.end(), PiZeroPtSorter());
00063   std::vector<const RecoTauPiZero*>::const_iterator wimp =
00064       signalPiZeros.begin();
00065 
00066   // Loop until we have a sturdy enough pizero
00067   reco::Candidate::LorentzVector totalP4 = tau.p4();
00068   while(wimp != signalPiZeros.end() && filter(*wimp, totalP4)) {
00069     totalP4 -= (*wimp)->p4();
00070     ++wimp;
00071   }
00072 
00073   if (wimp != signalPiZeros.begin()) {
00074     // We filtered stuff, update our tau
00075     double ptDifference = (totalP4 - tau.p4()).pt();
00076     tau.setisolationPFGammaCandsEtSum(
00077         tau.isolationPFGammaCandsEtSum() + ptDifference);
00078 
00079     // Update four vector
00080     tau.setP4(totalP4);
00081 
00082     std::vector<RecoTauPiZero> toMove;
00083     std::vector<RecoTauPiZero> newSignal;
00084     std::vector<RecoTauPiZero> newIsolation;
00085 
00086     // Build our new objects
00087     for (std::vector<const RecoTauPiZero*>::const_iterator iter =
00088          signalPiZeros.begin(); iter != wimp; ++iter) {
00089       toMove.push_back(**iter);
00090     }
00091     for (std::vector<const RecoTauPiZero*>::const_iterator iter =
00092          wimp; iter != signalPiZeros.end(); ++iter) {
00093       newSignal.push_back(**iter);
00094     }
00095     // Build our new isolation collection
00096     std::copy(toMove.begin(), toMove.end(), std::back_inserter(newIsolation));
00097     std::copy(tau.isolationPiZeroCandidates().begin(),
00098               tau.isolationPiZeroCandidates().end(),
00099               std::back_inserter(newIsolation));
00100 
00101     // Set the collections in the taus.
00102     tau.setsignalPiZeroCandidates(newSignal);
00103     tau.setisolationPiZeroCandidates(newIsolation);
00104 
00105     // Now we need to deal with the gamma candidates underlying moved pizeros.
00106     std::vector<PFCandidatePtr> pfcandsToMove = flattenPiZeros(toMove);
00107 
00108     // Copy the keys to move
00109     std::set<size_t> keysToMove;
00110     BOOST_FOREACH(const PFCandidatePtr& ptr, pfcandsToMove) {
00111       keysToMove.insert(ptr.key());
00112     }
00113 
00114     PFCandidateRefVector newSignalPFGammas;
00115     PFCandidateRefVector newSignalPFCands;
00116     PFCandidateRefVector newIsolationPFGammas = tau.isolationPFGammaCands();
00117     PFCandidateRefVector newIsolationPFCands = tau.isolationPFCands();
00118 
00119     // Move the necessary signal pizeros - what a mess!
00120     BOOST_FOREACH(const PFCandidateRef& ref, tau.signalPFCands()) {
00121       if (keysToMove.count(ref.key()))
00122         newIsolationPFCands.push_back(ref);
00123       else
00124         newSignalPFCands.push_back(ref);
00125     }
00126 
00127     BOOST_FOREACH(const PFCandidateRef& ref, tau.signalPFGammaCands()) {
00128       if (keysToMove.count(ref.key()))
00129         newIsolationPFGammas.push_back(ref);
00130       else
00131         newSignalPFGammas.push_back(ref);
00132     }
00133 
00134     tau.setsignalPFCands(newSignalPFCands);
00135     tau.setsignalPFCands(newSignalPFGammas);
00136     tau.setisolationPFGammaCands(newIsolationPFGammas);
00137     tau.setisolationPFCands(newIsolationPFCands);
00138   }
00139 }
00140 }}  // end namespace reco::tau
00141 #include "FWCore/Framework/interface/MakerMacros.h"
00142 DEFINE_EDM_PLUGIN(RecoTauModifierPluginFactory,
00143     reco::tau::RecoTauPhotonFilter,
00144     "RecoTauPhotonFilter");