CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoTauTag/RecoTau/plugins/RecoTauPiZeroTrivialPlugin.cc

Go to the documentation of this file.
00001 /*
00002  * RecoTauPiZeroTrivialPlugin
00003  *
00004  * Author: Evan K. Friis (UC Davis)
00005  *
00006  * Given an input PFJet, produces collection of trivial 'un-merged' PiZero
00007  * RecoTauPiZeros.  Each PiZero is composed of only one photon from
00008  * the jet.
00009  *
00010  * $Id $
00011  *
00012  */
00013 
00014 #include "RecoTauTag/RecoTau/interface/RecoTauPiZeroPlugins.h"
00015 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00016 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00017 #include "DataFormats/JetReco/interface/PFJet.h"
00018 #include "DataFormats/TauReco/interface/RecoTauPiZero.h"
00019 
00020 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00021 #include "RecoTauTag/RecoTau/interface/RecoTauQualityCuts.h"
00022 
00023 #include <boost/foreach.hpp>
00024 
00025 namespace reco { namespace tau {
00026 
00027 class RecoTauPiZeroTrivialPlugin : public RecoTauPiZeroBuilderPlugin {
00028   public:
00029     explicit RecoTauPiZeroTrivialPlugin(const edm::ParameterSet& pset);
00030     ~RecoTauPiZeroTrivialPlugin() {}
00031     return_type operator()(const reco::PFJet& jet) const;
00032   private:
00033     RecoTauQualityCuts qcuts_;
00034 };
00035 
00036 RecoTauPiZeroTrivialPlugin::RecoTauPiZeroTrivialPlugin(
00037     const edm::ParameterSet& pset):RecoTauPiZeroBuilderPlugin(pset),
00038     qcuts_(pset.getParameter<edm::ParameterSet>("qualityCuts")) {}
00039 
00040 RecoTauPiZeroBuilderPlugin::return_type RecoTauPiZeroTrivialPlugin::operator()(
00041     const reco::PFJet& jet) const {
00042   std::vector<PFCandidatePtr> pfGammaCands = qcuts_.filterRefs(pfGammas(jet));
00043   PiZeroVector output;
00044   output.reserve(pfGammaCands.size());
00045 
00046   BOOST_FOREACH(const PFCandidatePtr& gamma, pfGammaCands) {
00047     std::auto_ptr<RecoTauPiZero> piZero(new RecoTauPiZero(
00048             0, (*gamma).p4(), (*gamma).vertex(), 22, 1000, true,
00049             RecoTauPiZero::kTrivial));
00050     piZero->addDaughter(gamma);
00051     output.push_back(piZero);
00052   }
00053   return output.release();
00054 }
00055 
00056 }} // end reco::tau
00057 
00058 #include "FWCore/Framework/interface/MakerMacros.h"
00059 DEFINE_EDM_PLUGIN(RecoTauPiZeroBuilderPluginFactory,
00060     reco::tau::RecoTauPiZeroTrivialPlugin,
00061     "RecoTauPiZeroTrivialPlugin");