CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoTauTag/RecoTau/plugins/RecoTauCleaner.cc

Go to the documentation of this file.
00001 /*
00002  * RecoTauCleaner
00003  *
00004  * Author: Evan K. Friis, UC Davis
00005  *
00006  * Given a input collection of PFTaus, produces an output collection of PFTaus
00007  * such that no two PFTaus come from the same PFJet.  If multiple taus in the
00008  * collection come from the same PFJet, (dirty) they are ordered according to a
00009  * list of cleaners.  Each cleaner is a RecoTauCleanerPlugin, and returns a
00010  * double corresponding to the 'quality' for a given tau - an example would be
00011  * the level of isolation.  The set of dirty taus is then ranked
00012  * lexicographically by these cleaners, and the best one is placed in the
00013  * output collection.
00014  *
00015  * $Id $
00016  */
00017 
00018 #include <boost/ptr_container/ptr_vector.hpp>
00019 #include <algorithm>
00020 #include <memory>
00021 
00022 #include "FWCore/Framework/interface/EDProducer.h"
00023 #include "FWCore/Framework/interface/EventSetup.h"
00024 #include "FWCore/Framework/interface/ESHandle.h"
00025 #include "FWCore/Framework/interface/Event.h"
00026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00027 
00028 #include "RecoTauTag/RecoTau/interface/RecoTauBuilderPlugins.h"
00029 #include "RecoTauTag/RecoTau/interface/RecoTauCleaningTools.h"
00030 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00031 
00032 #include "DataFormats/TauReco/interface/PFTau.h"
00033 #include "DataFormats/TauReco/interface/PFTauFwd.h"
00034 
00035 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
00036 
00037 template<typename Prod>
00038 class RecoTauCleanerImpl : public edm::EDProducer {
00039   typedef reco::tau::RecoTauCleanerPlugin Cleaner;
00040   typedef boost::ptr_vector<Cleaner> CleanerList;
00041   // Define our output type - i.e. reco::PFTau OR reco::PFTauRef
00042   typedef typename Prod::value_type output_type;
00043 
00044   // Predicate that determines if two taus 'overlap' i.e. share a base PFJet
00045   class RemoveDuplicateJets {
00046     public:
00047       bool operator()(const reco::PFTauRef& a, const reco::PFTauRef& b)
00048       { return a->jetRef() == b->jetRef(); }
00049   };
00050 
00051   public:
00052     explicit RecoTauCleanerImpl(const edm::ParameterSet& pset);
00053     ~RecoTauCleanerImpl() {}
00054     void produce(edm::Event& evt, const edm::EventSetup& es);
00055 
00056   private:
00057     edm::InputTag tauSrc_;
00058     CleanerList cleaners_;
00059     // Optional selection on the output of the taus
00060     std::auto_ptr<StringCutObjectSelector<reco::PFTau> > outputSelector_;
00061 };
00062 
00063 template<typename Prod>
00064 RecoTauCleanerImpl<Prod>::RecoTauCleanerImpl(const edm::ParameterSet& pset) {
00065   tauSrc_ = pset.getParameter<edm::InputTag>("src");
00066   // Build our list of quality plugins
00067   typedef std::vector<edm::ParameterSet> VPSet;
00068   // Get each of our tau builders
00069   const VPSet& cleaners = pset.getParameter<VPSet>("cleaners");
00070   for (VPSet::const_iterator cleanerPSet = cleaners.begin();
00071       cleanerPSet != cleaners.end(); ++cleanerPSet) {
00072     // Get plugin name
00073     const std::string& pluginType =
00074       cleanerPSet->getParameter<std::string>("plugin");
00075     // Build the plugin
00076     cleaners_.push_back(
00077         RecoTauCleanerPluginFactory::get()->create(pluginType, *cleanerPSet));
00078   }
00079 
00080   // Check if we want to apply a final output selection
00081   if (pset.exists("outputSelection")) {
00082     std::string selection = pset.getParameter<std::string>("outputSelection");
00083     if (selection != "") {
00084       outputSelector_.reset(
00085           new StringCutObjectSelector<reco::PFTau>(selection));
00086     }
00087   }
00088 
00089   // Build the predicate that ranks our taus.  
00090   produces<Prod>();
00091 }
00092 
00093 namespace {
00094 // Template to convert a ref to desired output type
00095 template<typename T> const T convert(const reco::PFTauRef &tau);
00096 
00097 template<> const edm::RefToBase<reco::PFTau>
00098 convert<edm::RefToBase<reco::PFTau> >(const reco::PFTauRef &tau) {
00099   return edm::RefToBase<reco::PFTau>(tau);
00100 }
00101 
00102 template<> const reco::PFTauRef
00103 convert<reco::PFTauRef>(const reco::PFTauRef &tau) { return tau; }
00104 
00105 template<> const reco::PFTau
00106 convert<reco::PFTau>(const reco::PFTauRef &tau) { return *tau; }
00107 }
00108 
00109 template<typename Prod>
00110 void RecoTauCleanerImpl<Prod>::produce(edm::Event& evt,
00111                                    const edm::EventSetup& es) {
00112   // Update all our cleaners with the event info if they need it
00113   for (CleanerList::iterator cleaner = cleaners_.begin();
00114       cleaner != cleaners_.end(); ++cleaner) {
00115     cleaner->setup(evt, es);
00116   }
00117 
00118   // Get the input collection to clean
00119   edm::Handle<reco::CandidateView> input;
00120   evt.getByLabel(tauSrc_, input);
00121 
00122   // Cast the input candidates to Refs to real taus
00123   reco::PFTauRefVector inputRefs =
00124       reco::tau::castView<reco::PFTauRefVector>(input);
00125 
00126   // Make an STL algorithm friendly vector of refs
00127   typedef std::vector<reco::PFTauRef> PFTauRefs;
00128   // Collection of all taus. Some are from the same PFJet. We must clean them.
00129   PFTauRefs dirty(inputRefs.size());
00130  
00131   // Sort the input tau refs according to our predicate
00132   auto N = inputRefs.size();
00133   auto CN = cleaners_.size();
00134   int index[N];
00135   for (decltype(N) i=0; i!=N; ++i) index[i]=i;
00136   float  ranks[N*CN];
00137   // auto ranks = [rr,CN](int i, int j){return rr[i*CN+j];};
00138   // float const * *  rr = ranks;
00139   decltype(N) ii=0;
00140   for (decltype(N) ir=0; ir!=N; ++ir)
00141     for(decltype(N) cl=0;cl!=CN;++cl)
00142       ranks[ii++]=cleaners_[cl](inputRefs[ir]);
00143   assert(ii==(N*CN));
00144   const float * rr = ranks;
00145   std::sort(index,index+N,[rr,CN](int i, int j) { return std::lexicographical_compare(rr+i*CN,rr+i*CN+CN,rr+j*CN,rr+j*CN+CN);});
00146   for (decltype(N) ir=0; ir!=N; ++ir)
00147     dirty[ir]=inputRefs[index[ir]];
00148 
00149   // Clean the taus, ensuring that only one tau per jet is produced
00150   PFTauRefs cleanTaus = reco::tau::cleanOverlaps<PFTauRefs,
00151             RemoveDuplicateJets>(dirty);
00152 
00153   // create output collection
00154   std::auto_ptr<Prod> output(new Prod());
00155   //output->reserve(cleanTaus.size());
00156 
00157   // Copy clean refs into output
00158   for (PFTauRefs::const_iterator tau = cleanTaus.begin();
00159        tau != cleanTaus.end(); ++tau) {
00160     // If we are applying an output selection, check if it passes
00161     bool selected = true;
00162     if (outputSelector_.get() && !(*outputSelector_)(**tau)) {
00163       selected = false;
00164     }
00165     if (selected) {
00166       output->push_back(convert<output_type>(*tau));
00167     }
00168   }
00169   evt.put(output);
00170 }
00171 
00172 typedef RecoTauCleanerImpl<reco::PFTauCollection> RecoTauCleaner;
00173 typedef RecoTauCleanerImpl<reco::PFTauRefVector> RecoTauRefCleaner;
00174 
00175 #include "FWCore/Framework/interface/MakerMacros.h"
00176 DEFINE_FWK_MODULE(RecoTauCleaner);
00177 DEFINE_FWK_MODULE(RecoTauRefCleaner);
00178