CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTauTag/TauTagTools/plugins/PFTauViewRefSelector.cc

Go to the documentation of this file.
00001 /* \class CandViewRefSelector
00002  *
00003  * PFTau Selector based on a configurable cut.
00004  * Reads a edm::View<PFTau> as input
00005  * and saves a vector of references
00006  * Usage:
00007  *
00008  * module selectedCands = PFTauViewRefSelector {
00009  *   InputTag src = myCollection
00010  *   string cut = "pt > 15.0"
00011  * };
00012  *
00013  * \author: Luca Lista, INFN
00014  * \modifications: Evan Friis, UC Davis
00015  *
00016  */
00017 
00018 #include <boost/foreach.hpp>
00019 
00020 #include "FWCore/Framework/interface/MakerMacros.h"
00021 #include "CommonTools/UtilAlgos/interface/StringCutObjectSelector.h"
00022 
00023 #include "DataFormats/TauReco/interface/PFTauFwd.h"
00024 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00025 
00026 #include "FWCore/Framework/interface/EDFilter.h"
00027 #include "FWCore/Framework/interface/EventSetup.h"
00028 #include "FWCore/Framework/interface/ESHandle.h"
00029 #include "FWCore/Framework/interface/Event.h"
00030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00031 
00032 class PFTauViewRefSelector : public edm::EDFilter {
00033   public:
00034     explicit PFTauViewRefSelector(const edm::ParameterSet &pset);
00035     ~PFTauViewRefSelector() {}
00036     bool filter(edm::Event& evt, const edm::EventSetup& es);
00037   private:
00038     edm::InputTag src_;
00039     std::string cut_;
00040     std::auto_ptr<StringCutObjectSelector<reco::PFTau> > outputSelector_;
00041     bool filter_;
00042 };
00043 
00044 PFTauViewRefSelector::PFTauViewRefSelector(const edm::ParameterSet &pset) {
00045   src_ = pset.getParameter<edm::InputTag>("src");
00046   std::string cut = pset.getParameter<std::string>("cut");
00047   filter_ = pset.exists("filter") ? pset.getParameter<bool>("filter") : false;
00048   outputSelector_.reset(new StringCutObjectSelector<reco::PFTau>(cut));
00049   produces<reco::PFTauRefVector>();
00050 }
00051 
00052 bool
00053 PFTauViewRefSelector::filter(edm::Event& evt, const edm::EventSetup& es) {
00054   std::auto_ptr<reco::PFTauRefVector> output(new reco::PFTauRefVector());
00055   // Get the input collection to clean
00056   edm::Handle<reco::CandidateView> input;
00057   evt.getByLabel(src_, input);
00058   // Cast the input candidates to Refs to real taus
00059   reco::PFTauRefVector inputRefs =
00060       reco::tau::castView<reco::PFTauRefVector>(input);
00061 
00062   BOOST_FOREACH(const reco::PFTauRef &tau, inputRefs) {
00063     if (outputSelector_.get() && (*outputSelector_)(*tau)) {
00064       output->push_back(tau);
00065     }
00066   }
00067   size_t outputSize = output->size();
00068   evt.put(output);
00069   // Filter if desired and no objects passed our cut
00070   return !(filter_ && outputSize == 0);
00071 }
00072 
00073 DEFINE_FWK_MODULE(PFTauViewRefSelector);
00074 
00075 #include "CommonTools/UtilAlgos/interface/SingleObjectSelector.h"
00076 #include "DataFormats/JetReco/interface/PFJet.h"
00077 typedef SingleObjectSelector<
00078           edm::View<reco::Jet>,
00079           StringCutObjectSelector<reco::Jet, true>
00080        > JetViewRefSelector;
00081 
00082 DEFINE_FWK_MODULE(JetViewRefSelector);