CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoTauTag/TauTagTools/plugins/RecoTauDistanceFromTruthPlugin.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/RecoTau/interface/RecoTauBuilderPlugins.h"
00002 
00003 #include "DataFormats/TauReco/interface/PFTau.h"
00004 #include "DataFormats/TauReco/interface/PFTauFwd.h"
00005 #include "DataFormats/JetReco/interface/GenJet.h"
00006 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00007 #include "DataFormats/Common/interface/Association.h"
00008 //#include "DataFormats/Common/interface/AssociativeIterator.h"
00009 
00010 #include <boost/foreach.hpp>
00011 #include <boost/bind.hpp>
00012 #include <boost/iterator/filter_iterator.hpp>
00013 
00014 namespace tautools {
00015 
00016 class RecoTauDistanceFromTruthPlugin : public reco::tau::RecoTauCleanerPlugin {
00017   public:
00018     RecoTauDistanceFromTruthPlugin(const edm::ParameterSet& pset);
00019     virtual ~RecoTauDistanceFromTruthPlugin() {}
00020     double operator()(const reco::PFTauRef&) const;
00021     void beginEvent();
00022   private:
00023     edm::InputTag matchingSrc_;
00024     typedef edm::Association<reco::GenJetCollection> GenJetAssociation;
00025     edm::Handle<GenJetAssociation> genTauMatch_;
00026 };
00027 
00028 RecoTauDistanceFromTruthPlugin::RecoTauDistanceFromTruthPlugin(
00029     const edm::ParameterSet& pset): reco::tau::RecoTauCleanerPlugin(pset) {
00030   matchingSrc_ = pset.getParameter<edm::InputTag>("matching");
00031 }
00032 
00033 void RecoTauDistanceFromTruthPlugin::beginEvent() {
00034   // Load the matching information
00035   evt()->getByLabel(matchingSrc_, genTauMatch_);
00036 }
00037 
00038 // Helpers
00039 namespace {
00040   // Returns the squared momentum difference between two candidates
00041   double momentumDifference(const reco::Candidate* candA,
00042       const reco::Candidate* candB) {
00043     reco::Candidate::LorentzVector difference =
00044       candA->p4() - candB->p4();
00045     return difference.P2();
00046   }
00047 
00048   // Finds the best match for an input <cand> from an input colleciton.
00049   // Only objects with matching charge are considered.  The best match
00050   // has the lowest [momentumDifference] with the input <cand>
00051   template<typename InputIterator>
00052   InputIterator findBestMatch(const reco::Candidate* cand,
00053       InputIterator begin, InputIterator end) {
00054 
00055     typedef const reco::Candidate* CandPtr;
00056     using boost::bind;
00057     using boost::function;
00058     using boost::filter_iterator;
00059     // Build a charge matching function
00060     typedef function<bool (CandPtr)> CandPtrBoolFn;
00061     CandPtrBoolFn chargeMatcher =
00062       bind(&reco::Candidate::charge, cand) == bind(&reco::Candidate::charge, _1);
00063 
00064     // Only match those objects that have the same charge
00065     typedef filter_iterator<CandPtrBoolFn, InputIterator> Iterator;
00066     Iterator begin_filtered(chargeMatcher, begin, end);
00067     Iterator end_filtered(chargeMatcher, end, end);
00068 
00069     Iterator result = std::min_element(begin_filtered, end_filtered,
00070         momentumDifference);
00071     return result.base();
00072   }
00073 } // end anon. namespace
00074 
00075 double RecoTauDistanceFromTruthPlugin::operator()(const reco::PFTauRef& tauRef) const {
00076 
00077   GenJetAssociation::reference_type truth = (*genTauMatch_)[tauRef];
00078 
00079   // Check if the matching exists, if not return +infinity
00080   if (truth.isNull())
00081     return std::numeric_limits<double>::infinity();
00082 
00083   // screw all this noise
00084   return std::abs(tauRef->pt() - truth->pt());
00085 
00086   typedef const reco::Candidate* CandPtr;
00087   typedef std::set<CandPtr> CandPtrSet;
00088   typedef std::vector<CandPtr> CandPtrVector;
00089   typedef std::list<CandPtr> CandPtrList;
00090   // Store all of our reco and truth candidates
00091   CandPtrList recoCands;
00092   CandPtrSet truthCandSet;
00093 
00094   BOOST_FOREACH(const reco::RecoTauPiZero& piZero,
00095       tauRef->signalPiZeroCandidates()) {
00096     recoCands.push_back(&piZero);
00097   }
00098 
00099   BOOST_FOREACH(const reco::PFCandidateRef& pfch,
00100       tauRef->signalPFChargedHadrCands()) {
00101     recoCands.push_back(pfch.get());
00102   }
00103 
00104   // Use a set to contain the truth pointers so we ensure that no pizeros
00105   // are entered twice.
00106   BOOST_FOREACH(const reco::CandidatePtr& ptr,
00107       truth->daughterPtrVector()) {
00108     // Get mother pi zeros in the case of gammas
00109     if (ptr->pdgId() == 22)
00110       truthCandSet.insert(ptr->mother());
00111     else
00112       truthCandSet.insert(ptr.get());
00113   }
00114 
00115   //Convert truth cands from set to vector so we can sort it.
00116   CandPtrVector truthCands(truthCandSet.begin(), truthCandSet.end());
00117 
00118   // Sort the truth candidates by descending pt
00119   std::sort(truthCands.begin(), truthCands.end(),
00120       boost::bind(&reco::Candidate::pt, _1) > boost::bind(&reco::Candidate::pt, _2));
00121 
00122   double quality = 0.0;
00123   BOOST_FOREACH(CandPtr truthCand, truthCands) {
00124     // Find the best reco match for this truth cand
00125     CandPtrList::iterator recoMatch = findBestMatch(truthCand,
00126         recoCands.begin(), recoCands.end());
00127 
00128     // Check if this truth cand is matched
00129     if (recoMatch != recoCands.end()) {
00130       // Add a penalty factor based on how different the reconstructed
00131       // is w.r.t. the true momentum
00132       quality += momentumDifference(truthCand, *recoMatch);
00133       // Remove this reco cand from future matches
00134       recoCands.erase(recoMatch);
00135     } else {
00136       // this truth cand was not matched!
00137       quality += truthCand->p4().P2();
00138     }
00139   }
00140 
00141   // Now add a penalty for the unmatched reco stuff
00142   BOOST_FOREACH(CandPtr recoCand, recoCands) {
00143     quality += recoCand->p4().P2();
00144   }
00145 
00146   return quality;
00147 }
00148 
00149 } // end tautools namespace
00150 
00151 
00152 // Register our plugin
00153 #include "FWCore/Framework/interface/MakerMacros.h"
00154 DEFINE_EDM_PLUGIN(RecoTauCleanerPluginFactory, tautools::RecoTauDistanceFromTruthPlugin, "RecoTauDistanceFromTruthPlugin");