CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/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.getParameterSet(
00039           "qualityCuts").getParameterSet("signalQualityCuts")) {}
00040 
00041 RecoTauPiZeroBuilderPlugin::return_type RecoTauPiZeroTrivialPlugin::operator()(
00042     const reco::PFJet& jet) const {
00043   std::vector<PFCandidatePtr> pfGammaCands = qcuts_.filterRefs(pfGammas(jet));
00044   PiZeroVector output;
00045   output.reserve(pfGammaCands.size());
00046 
00047   BOOST_FOREACH(const PFCandidatePtr& gamma, pfGammaCands) {
00048     std::auto_ptr<RecoTauPiZero> piZero(new RecoTauPiZero(
00049             0, (*gamma).p4(), (*gamma).vertex(), 22, 1000, true,
00050             RecoTauPiZero::kTrivial));
00051     piZero->addDaughter(gamma);
00052     output.push_back(piZero);
00053   }
00054   return output.release();
00055 }
00056 
00057 }} // end reco::tau
00058 
00059 #include "FWCore/Framework/interface/MakerMacros.h"
00060 DEFINE_EDM_PLUGIN(RecoTauPiZeroBuilderPluginFactory,
00061     reco::tau::RecoTauPiZeroTrivialPlugin,
00062     "RecoTauPiZeroTrivialPlugin");