00001 /* 00002 * RecoTauPiZeroStripPlugin 00003 * 00004 * Merges PFGammas in a PFJet into Candidate piZeros defined as 00005 * strips in eta-phi. 00006 * 00007 * Author: Michail Bachtis (University of Wisconsin) 00008 * 00009 * Code modifications: Evan Friis (UC Davis) 00010 * 00011 * $Id $ 00012 */ 00013 #include <algorithm> 00014 #include <memory> 00015 00016 #include "RecoTauTag/RecoTau/interface/RecoTauPiZeroPlugins.h" 00017 00018 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h" 00019 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h" 00020 #include "DataFormats/VertexReco/interface/Vertex.h" 00021 #include "DataFormats/VertexReco/interface/VertexFwd.h" 00022 #include "DataFormats/TauReco/interface/RecoTauPiZero.h" 00023 #include "DataFormats/JetReco/interface/PFJet.h" 00024 #include "CommonTools/CandUtils/interface/AddFourMomenta.h" 00025 #include "DataFormats/Math/interface/deltaPhi.h" 00026 00027 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h" 00028 #include "RecoTauTag/RecoTau/interface/RecoTauQualityCuts.h" 00029 00030 namespace reco { namespace tau { 00031 00032 class RecoTauPiZeroStripPlugin : public RecoTauPiZeroBuilderPlugin { 00033 public: 00034 explicit RecoTauPiZeroStripPlugin(const edm::ParameterSet& pset); 00035 virtual ~RecoTauPiZeroStripPlugin() {} 00036 // Return type is auto_ptr<PiZeroVector> 00037 return_type operator()(const reco::PFJet& jet) const; 00038 // Hook to update PV information 00039 virtual void beginEvent(); 00040 00041 private: 00042 // PV needed for quality cuts 00043 edm::InputTag pvSrc_; 00044 RecoTauQualityCuts qcuts_; 00045 00046 std::vector<int> inputPdgIds_; //type of candidates to clusterize 00047 double etaAssociationDistance_;//eta Clustering Association Distance 00048 double phiAssociationDistance_;//phi Clustering Association Distance 00049 00050 AddFourMomenta p4Builder_; 00051 }; 00052 00053 RecoTauPiZeroStripPlugin::RecoTauPiZeroStripPlugin( 00054 const edm::ParameterSet& pset):RecoTauPiZeroBuilderPlugin(pset), 00055 qcuts_(pset.getParameter<edm::ParameterSet>("qualityCuts")) 00056 { 00057 pvSrc_ = pset.getParameter<edm::InputTag>("primaryVertexSrc"); 00058 inputPdgIds_ = pset.getParameter<std::vector<int> >( 00059 "stripCandidatesParticleIds"); 00060 etaAssociationDistance_ = pset.getParameter<double>( 00061 "stripEtaAssociationDistance"); 00062 phiAssociationDistance_ = pset.getParameter<double>( 00063 "stripPhiAssociationDistance"); 00064 } 00065 00066 // Update the primary vertex 00067 void RecoTauPiZeroStripPlugin::beginEvent() { 00068 edm::Handle<reco::VertexCollection> pvHandle; 00069 evt()->getByLabel(pvSrc_, pvHandle); 00070 if (pvHandle->size()) { 00071 qcuts_.setPV(reco::VertexRef(pvHandle, 0)); 00072 } 00073 } 00074 00075 RecoTauPiZeroStripPlugin::return_type RecoTauPiZeroStripPlugin::operator()( 00076 const reco::PFJet& jet) const { 00077 // Get list of gamma candidates 00078 typedef std::vector<reco::PFCandidatePtr> PFCandPtrs; 00079 typedef PFCandPtrs::iterator PFCandIter; 00080 PiZeroVector output; 00081 00082 // Get the candidates passing our quality cuts 00083 //PFCandPtrs candsVector = qcuts_.filterRefs(pfCandidates(jet, inputPdgIds_)); 00084 PFCandPtrs candsVector = qcuts_.filterRefs(pfGammas(jet)); 00085 00086 // Convert to stl::list to allow fast deletions 00087 typedef std::list<reco::PFCandidatePtr> PFCandPtrList; 00088 typedef std::list<reco::PFCandidatePtr>::iterator PFCandPtrListIter; 00089 PFCandPtrList cands; 00090 cands.insert(cands.end(), candsVector.begin(), candsVector.end()); 00091 00092 while (cands.size() > 0) { 00093 // Seed this new strip, and delete it from future strips 00094 PFCandidatePtr seed = cands.front(); 00095 cands.pop_front(); 00096 00097 // Add a new candidate to our collection using this seed 00098 std::auto_ptr<RecoTauPiZero> strip(new RecoTauPiZero( 00099 *seed, RecoTauPiZero::kStrips)); 00100 strip->addDaughter(seed); 00101 00102 // Find all other objects in the strip 00103 PFCandPtrListIter stripCand = cands.begin(); 00104 while(stripCand != cands.end()) { 00105 if( fabs(strip->eta() - (*stripCand)->eta()) < etaAssociationDistance_ 00106 && fabs(deltaPhi(*strip, **stripCand)) < phiAssociationDistance_ ) { 00107 // Add candidate to strip 00108 strip->addDaughter(*stripCand); 00109 // Update the strips four momenta 00110 p4Builder_.set(*strip); 00111 // Delete this candidate from future strips and move on to 00112 // the next potential candidate 00113 stripCand = cands.erase(stripCand); 00114 } else { 00115 // This candidate isn't compatabile - just move to the next candidate 00116 ++stripCand; 00117 } 00118 } 00119 output.push_back(strip); 00120 } 00121 return output.release(); 00122 } 00123 }} // end namespace reco::tau 00124 00125 #include "FWCore/Framework/interface/MakerMacros.h" 00126 DEFINE_EDM_PLUGIN(RecoTauPiZeroBuilderPluginFactory, 00127 reco::tau::RecoTauPiZeroStripPlugin, "RecoTauPiZeroStripPlugin");